source: trunk/Roms_agrif/zoombc_2D_patrick.F @ 4

Last change on this file since 4 was 4, checked in by pinsard, 17 years ago

modification according to /usr/temp/vech/quiberon/Roms_tools/Roms_Agrif_pisces/

File size: 38.3 KB
Line 
1!====================================================================
2!                   subroutine Agrif_u2dbc_interp_tile
3!====================================================================
4!   
5#include "cppdefs.h"     
6#ifdef AGRIF     
7      subroutine u2dbc_interp_tile(Istr,Iend,Jstr,Jend
8     &        ,DU_west4,DU_east4,DU_west6,DU_east6
9     &        ,DU_south4,DU_north4,DU_south6,DU_north6)
10      use AGRIF_Util
11!     implicit none
12# include "param.h"
13# include "grid.h"
14# include "ocean2d.h"   
15# include "scalars.h"
16# include "climat.h"
17# include "boundary.h"
18# include "zoom.h"
19# include "coupling.h"
20      integer Istr,Iend,Jstr,Jend, i,j
21      real t1,t2,t3,t7,t8,t9,t10,t11,t12,t13,t14,t15,t16,t17,t18,t19
22      logical ptinterp
23      real t4,t5,t6,tfin,c1,c2,c3,cff1,cff2,cff3
24      external :: u2dinterp
25      integer :: irhox, irhoy, irhot
26      real :: rrhox, rrhoy, rrhot
27      real :: tinterp, onemtinterp
28      integer :: iter
29      integer :: ipu, jpu, parentnnew
30      integer :: parentnbstep
31      real :: cffx, cffy
32      real,dimension(PRIVATE_1DETA_SCRATCH_ARRAY) :: DU_west, DU_east,
33     &  DU_west6,DU_east6
34      real,dimension(PRIVATE_1DETA_SCRATCH_ARRAY,0:NWEIGHT) :: 
35     &                      DU_west4, DU_east4
36      real,dimension(PRIVATE_1DXI_SCRATCH_ARRAY) :: DU_south,DU_north,
37     &                      DU_south6,DU_north6
38      real,dimension(PRIVATE_1DXI_SCRATCH_ARRAY,0:NWEIGHT) :: 
39     &                      DU_south4,DU_north4     
40      real :: lastcff
41!$AGRIF_DO_NOT_TREAT
42      real,dimension(:,:,:),pointer :: coarsevalues
43      integer :: nbgrid, indinterp
44      common/interp2d/coarsevalues, nbgrid, indinterp
45!$AGRIF_END_DO_NOT_TREAT     
46# ifdef MPI
47       include 'mpif.h'
48# endif
49 
50# ifdef MPI
51#  define LOCALLM Lmmpi
52#  define LOCALMM Mmmpi
53# else
54#  define LOCALLM Lm
55#  define LOCALMM Mm
56# endif   
57!         
58# include "compute_auxiliary_bounds.h"
59!     return
60      irhox=Agrif_Irhox()
61      irhoy=Agrif_Irhoy()
62      irhot=Agrif_Irhot()
63
64      rrhox = real(irhox)
65      rrhoy = real(irhoy)
66      rrhot = real(irhot)
67
68
69      parentnbstep=Agrif_Parent_Nb_Step()
70      if (U2DTimeindex .NE. parentnbstep) then
71
72        do J=JstrR,JendR
73        do i=IstrR,IendR
74        dinterp(i,j) = 0.
75        enddo
76        enddo
77             
78C$OMP BARRIER
79C$OMP MASTER     
80
81        tinterp=1.
82# ifdef MASKING
83        Agrif_UseSpecialValue = .true.
84# endif
85        Agrif_SpecialValue = 0.
86
87        nbgrid = Agrif_Fixed()
88C        indinterp = 0 indinterp is set to zero in zetabc
89           
90        Call Agrif_Bc_variable(dinterp,DU_avg2,calledweight=tinterp,
91     &    procname = u2dinterp)
92                 
93        Agrif_UseSpecialvalue=.false. 
94C$OMP END MASTER
95C$OMP BARRIER
96
97          lastcff=0.
98          do iter=iif,iif+irhot-1
99          lastcff=lastcff+((iter+1-iif)/rrhot)*
100     &        weight2(iif+irhot-1,iter)
101          enddo
102         
103          lastcff=1./lastcff
104               
105# ifdef AGRIF_OBC_WEST
106      if (WESTERN_EDGE) then
107          i = Istr
108
109          do j=Jstr,Jend
110          DU_west4(j,iif-1) = DU_west3(i,j,iif-1)
111          enddo
112         
113          if (iif+irhot-1<=nfast) then
114           DU_west6=0.
115          do iter=0,iif-1
116          cff1=weight2(iif+irhot-1,iter)
117            do j=Jstr,Jend
118             DU_west6(j) = DU_west6(j) +cff1*DU_west4(j,iter)
119            enddo
120          enddo
121         
122          do iter=iif,iif+irhot-2
123          do j=Jstr,Jend
124          DU_west6(j)=DU_west6(j)+((iif+irhot-1-iter)/rrhot)*
125     &        weight2(iif+irhot-1,iter)*DU_west4(j,iif-1)
126          enddo
127          enddo
128         
129
130         
131          do j=Jstr,Jend
132          DU_west6(j)=lastcff*((dinterp(i,j)/rrhoy)-DU_west6(j))
133#  ifdef MASKING
134     &                               * umask(i,j)
135#  endif                     
136          enddo     
137         
138          do j=Jstr,Jend
139          DU_west6(j)=(DU_west6(j)-DU_west4(j,iif-1))/rrhot
140          enddo 
141          endif
142                                       
143      endif
144# endif 
145# ifdef AGRIF_OBC_EAST
146      if (EASTERN_EDGE) then
147          i=Iend+1         
148         
149          do j=Jstr,Jend
150          DU_east4(j,iif-1) = DU_east3(i,j,iif-1)
151          enddo
152         
153          if (iif+irhot-1<=nfast) then
154           DU_east6=0.
155          do iter=0,iif-1
156          cff1=weight2(iif+irhot-1,iter)
157            do j=Jstr,Jend
158             DU_east6(j) = DU_east6(j) +cff1*DU_east4(j,iter)
159            enddo
160          enddo
161         
162          do iter=iif,iif+irhot-2
163          do j=Jstr,Jend
164          DU_east6(j)=DU_east6(j)+((iif+irhot-1-iter)/rrhot)*
165     &        weight2(iif+irhot-1,iter)*DU_east4(j,iif-1)
166          enddo
167          enddo
168
169          do j=Jstr,Jend
170          DU_east6(j)=lastcff*((dinterp(i,j)/rrhoy)-DU_east6(j))
171#  ifdef MASKING
172     &                               * umask(i,j)
173#  endif         
174          enddo     
175         
176          do j=Jstr,Jend
177          DU_east6(j)=(DU_east6(j)-DU_east4(j,iif-1))/rrhot
178          enddo 
179          endif         
180                                 
181      endif
182# endif 
183
184# ifdef AGRIF_OBC_SOUTH
185
186      if (SOUTHERN_EDGE) then
187          j=Jstr-1
188         
189          do i=Istr,Iend
190          DU_south4(i,iif-1) = DU_south3(i,j,iif-1)
191          enddo
192         
193          if (iif+irhot-1<=nfast) then
194           DU_south6=0.
195          do iter=0,iif-1
196          cff1=weight2(iif+irhot-1,iter)
197            do i=Istr,Iend
198             DU_south6(i) = DU_south6(i) +cff1*DU_south4(i,iter)
199            enddo
200          enddo
201         
202          do iter=iif,iif+irhot-2
203          do i=Istr,Iend
204          DU_south6(i)=DU_south6(i)+((iif+irhot-1-iter)/rrhot)*
205     &        weight2(iif+irhot-1,iter)*DU_south4(i,iif-1)
206          enddo
207          enddo
208         
209          do i=Istr,Iend
210          DU_south6(i)=lastcff*((dinterp(i,j)/rrhox)-DU_south6(i))
211#  ifdef MASKING
212     &                               * umask(i,j)
213#  endif         
214          enddo     
215         
216          do i=Istr,Iend
217          DU_south6(i)=(DU_south6(i)-DU_south4(i,iif-1))/rrhot
218          enddo 
219          endif           
220                                       
221      endif
222     
223# endif 
224# ifdef AGRIF_OBC_NORTH
225      if (NORTHERN_EDGE) then 
226          j=Jend+1
227         
228          do i=Istr,Iend
229          DU_north4(i,iif-1) = DU_north3(i,j,iif-1)
230          enddo
231         
232          if (iif+irhot-1<=nfast) then
233           DU_north6=0.
234          do iter=0,iif-1
235          cff1=weight2(iif+irhot-1,iter)
236            do i=Istr,Iend
237             DU_north6(i) = DU_north6(i) +cff1*DU_north4(i,iter)
238            enddo
239          enddo
240         
241          do iter=iif,iif+irhot-2
242          do i=Istr,Iend
243          DU_north6(i)=DU_north6(i)+((iif+irhot-1-iter)/rrhot)*
244     &        weight2(iif+irhot-1,iter)*DU_north4(i,iif-1)
245          enddo
246          enddo
247         
248
249         
250          do i=Istr,Iend
251          DU_north6(i)=lastcff*((dinterp(i,j)/rrhox)-DU_north6(i))
252#  ifdef MASKING
253     &                               * umask(i,j)
254#  endif             
255          enddo     
256         
257          do i=Istr,Iend
258          DU_north6(i)=(DU_north6(i)-DU_north4(i,iif-1))/rrhot
259          enddo 
260          endif
261
262      endif
263# endif
264        U2DTimeindex = parentnbstep
265      endif
266
267# ifdef AGRIF_OBC_WEST
268      if (WESTERN_EDGE) then
269          i = Istr       
270
271          do j=Jstr,Jend
272           DU_west4(j,iif-1) = DU_west3(i,j,iif-1)
273          enddo
274           
275          do j=Jstr,Jend
276            DU_west(j)=DU_west4(j,iif-1)+DU_west6(j)
277          enddo
278       endif
279# endif 
280
281# ifdef AGRIF_OBC_EAST
282      if (EASTERN_EDGE) then
283          i=Iend+1
284         
285          do j=Jstr,Jend
286           DU_east4(j,iif-1) = DU_east3(i,j,iif-1)
287          enddo
288           
289          do j=Jstr,Jend
290            DU_east(j)=DU_east4(j,iif-1)+DU_east6(j)
291          enddo 
292        endif
293#  endif
294
295# ifdef AGRIF_OBC_SOUTH
296      if (SOUTHERN_EDGE) then
297          j=Jstr-1
298         
299          do i=Istr,Iend
300           DU_south4(i,iif-1) = DU_south3(i,j,iif-1)
301          enddo
302           
303          do i=Istr,Iend
304            DU_south(i)=DU_south4(i,iif-1)+DU_south6(i)
305          enddo
306       endif
307# endif
308
309# ifdef AGRIF_OBC_NORTH
310      if (NORTHERN_EDGE) then
311          j=Jend+1
312
313          do i=Istr,Iend
314           DU_north4(i,iif-1) = DU_north3(i,j,iif-1)
315          enddo
316           
317          do i=Istr,Iend
318            DU_north(i)=DU_north4(i,iif-1)+DU_north6(i)
319          enddo         
320       endif
321# endif 
322         
323
324# if defined M2_FRC_BRY || defined M2NUDGING
325!
326! Apply the value to ubclm or ubarbry
327!
328      cffx = g*dtfast*2./(1.+rrhox)
329
330#  ifdef AGRIF_2WAY
331      cffx = 0.
332#  endif
333
334#  ifdef AGRIF_OBC_WEST
335      if (WESTERN_EDGE) then
336        do j=Jstr,Jend
337#   ifdef M2_FRC_BRY
338          ubarbry_west(j)=
339#   else
340          ubclm(Istr,j)= (cffx/om_u(Istr,j))*
341     &                   (SSH(Istr,j)-zeta(Istr,j,knew))  +
342#   endif
343#   ifdef AGRIF_FLUX_BC
344     &              (2.*DU_west(j)/((h(Istr-1,j)+zeta(Istr-1,j,knew)
345     &                                 +h(Istr,j)+zeta(Istr,j,knew))
346     &                                                *on_u(Istr,j)))
347#   else
348     &                  DU_west(j)
349#   endif
350#   ifdef MASKING
351     &                  *umask(Istr,j)
352#   endif
353        enddo   
354
355      endif
356#  endif
357
358#  ifdef AGRIF_OBC_EAST
359      if (EASTERN_EDGE) then
360        do j=Jstr,Jend
361#   ifdef M2_FRC_BRY
362          ubarbry_east(j)=
363#   else
364          ubclm(Iend+1,j)=-(cffx/om_u(Iend+1,j))*
365     &                     (SSH(Iend,j)-zeta(Iend,j,knew)) +
366#   endif
367#   ifdef AGRIF_FLUX_BC
368     &                  (2.*DU_east(j)/(( h(Iend,j)+zeta(Iend,j,knew)
369     &                              +h(Iend+1,j)+zeta(Iend+1,j,knew))
370     &                                              *on_u(Iend+1,j)))
371#   else
372     &                      DU_east(j)
373#   endif
374#   ifdef MASKING
375     &                      *umask(Iend+1,j)
376#   endif
377        enddo       
378      endif
379#  endif
380
381#  ifdef AGRIF_OBC_SOUTH
382      if (SOUTHERN_EDGE) then       
383        do i=Istr,Iend
384#   ifdef M2_FRC_BRY
385          ubarbry_south(i)=
386#   else
387          ubclm(i,Jstr-1)=
388#   endif
389#   ifdef AGRIF_FLUX_BC
390     &              (2.*DU_south(i)/(( h(i,Jstr-1)+zeta(i,Jstr-1,knew)
391     &                           +h(i-1,Jstr-1)+zeta(i-1,Jstr-1,knew))
392     &                                                *on_u(i,Jstr-1)))
393#   else
394     &                  DU_south(i)
395#   endif
396#   ifdef MASKING
397     &                  *umask(i,Jstr-1)
398#   endif
399      enddo
400      endif
401#  endif
402
403#  ifdef AGRIF_OBC_NORTH
404      if (NORTHERN_EDGE) then
405        do i=Istr,Iend
406#   ifdef M2_FRC_BRY
407          ubarbry_north(i)=
408#   else
409          ubclm(i,Jend+1)=
410#   endif
411#   ifdef AGRIF_FLUX_BC
412     &            (2.*DU_north(i)/(( h(i,Jend+1)+zeta(i,Jend+1,knew)
413     &                         +h(i-1,Jend+1)+zeta(i-1,Jend+1,knew))
414     &                                             *on_u(i,Jend+1)))
415#   else
416     &                DU_north(i)
417#   endif
418#   ifdef MASKING
419     &                *umask(i,Jend+1)
420#   endif
421        enddo
422
423      endif
424#  endif
425# endif /* M2_FRC_BRY || M2NUDGING */
426
427      return
428      end
429     
430      subroutine u2Dinterp(tabres,i1,i2,j1,j2)
431      implicit none
432# include "param.h"
433# include "grid.h"
434# include "ocean2d.h"
435# include "scalars.h"
436
437       integer i1,i2,j1,j2
438       real tabres(i1:i2,j1:j2)
439       integer i,j,iter, isize
440       real :: cff1
441!$AGRIF_DO_NOT_TREAT
442      real,dimension(:,:,:),pointer :: coarsevalues
443      integer :: nbgrid, indinterp
444      common/interp2d/coarsevalues, nbgrid, indinterp 
445!$AGRIF_END_DO_NOT_TREAT 
446
447       integer :: oldindinterp
448       
449       isize = (j2-j1+1)*(i2-i1+1)
450       
451       IF (iif == 1) THEN
452       IF (.NOT.Associated(coarsevalues)) THEN
453       Allocate(coarsevalues(isize,0:nfast,1))
454       ELSE       
455       CALL checksizeinterp(indinterp+isize,nfast)
456       ENDIF
457       
458       oldindinterp = indinterp
459       
460       do j=j1,j2
461       do i=max(i1,lbound(h,1)+1),i2
462       oldindinterp=oldindinterp+1
463       coarsevalues(oldindinterp,0,nbgrid) = 
464     &  0.5*(h(i-1,j)+zeta(i-1,j,kstp)+h(i,j)+
465     &  zeta(i,j,kstp))*on_u(i,j)*ubar(i,j,kstp)
466       enddo
467       enddo
468             
469       ENDIF
470       
471       oldindinterp = indinterp
472       
473       do j=j1,j2
474       do i=max(i1,lbound(h,1)+1),i2
475       oldindinterp=oldindinterp+1
476       coarsevalues(oldindinterp,iif,nbgrid) =
477     &   0.5*(h(i-1,j)+zeta(i-1,j,knew)+h(i,j)+
478     &  zeta(i,j,knew))*on_u(i,j)*ubar(i,j,knew)
479       enddo
480       enddo
481         
482       do iter=0,iif
483          cff1 = weight2(iif,iter)
484          call copy1d(tabres,coarsevalues(indinterp+1,iter,nbgrid),
485     &                cff1,isize)
486       enddo
487       
488       indinterp = oldindinterp     
489       
490      return
491      end     
492!
493!====================================================================
494!                   subroutine Agrif_v2dbc_interp_tile
495!====================================================================
496!
497      subroutine v2dbc_interp_tile(Istr,Iend,Jstr,Jend
498     &        ,DV_west4,DV_east4,DV_west6,DV_east6
499     &        ,DV_south4,DV_north4,DV_south6,DV_north6)
500      use AGRIF_Util
501!     implicit none
502# include "param.h"
503# include "grid.h"
504# include "ocean2d.h"
505# include "scalars.h"
506# include "climat.h"
507# include "boundary.h"
508# include "zoom.h"
509# include "coupling.h"
510      integer Istr,Iend,Jstr,Jend, i,j
511      real t1,t2,t3,t4,t5,t6,dv1,dv1np1,dv2,tfin,c1,c2,c3
512      real t7,t8,t9,t10,t11
513      external :: v2dinterp
514      integer :: irhox, irhoy, irhot
515      real :: rrhox, rrhoy, rrhot
516      real :: tinterp, onemtinterp
517      integer :: iter 
518      integer :: parentnbstep   
519      real :: cffy
520      real,dimension(PRIVATE_1DETA_SCRATCH_ARRAY) :: DV_west, DV_east,
521     &  DV_west6,DV_east6
522      real,dimension(PRIVATE_1DETA_SCRATCH_ARRAY,0:NWEIGHT) :: 
523     &                      DV_west4, DV_east4
524      real,dimension(PRIVATE_1DXI_SCRATCH_ARRAY) :: DV_south,DV_north,
525     &                      DV_south6,DV_north6
526      real,dimension(PRIVATE_1DXI_SCRATCH_ARRAY,0:NWEIGHT) :: 
527     &                      DV_south4,DV_north4 
528
529      real :: lastcff
530!$AGRIF_DO_NOT_TREAT
531      real,dimension(:,:,:),pointer :: coarsevalues
532      integer :: nbgrid, indinterp
533      common/interp2d/coarsevalues, nbgrid, indinterp
534!$AGRIF_END_DO_NOT_TREAT     
535#ifdef MPI
536       include 'mpif.h'
537#endif
538 
539# ifdef MPI
540#  define LOCALLM Lmmpi
541#  define LOCALMM Mmmpi
542# else
543#  define LOCALLM Lm
544#  define LOCALMM Mm
545# endif   
546!         
547# include "compute_auxiliary_bounds.h"
548!
549!      return
550      irhox=Agrif_Irhox()
551      irhoy=Agrif_Irhoy()
552      irhot=Agrif_Irhot()
553     
554      rrhox = real(irhox)
555      rrhoy = real(irhoy)
556      rrhot = real(irhot)
557
558
559      parentnbstep=Agrif_Parent_Nb_Step()
560     
561      if (V2DTimeindex .NE. parentnbstep) then
562
563       
564        do J=JstrR,JendR
565        do i=IstrR,IendR
566        dinterp(i,j) = 0.
567        enddo
568        enddo
569
570C$OMP BARRIER
571C$OMP MASTER
572
573        tinterp=1.
574#ifdef MASKING       
575        Agrif_UseSpecialValue = .true.
576#endif
577
578        Agrif_SpecialValue = 0.
579       
580        nbgrid = Agrif_Fixed()
581C        indinterp = 0 indinterp is set to zero in zetabc
582        Call Agrif_Bc_variable(dinterp,DV_avg2,calledweight=tinterp,
583     &   procname = v2dinterp)
584                 
585        Agrif_UseSpecialvalue=.false.
586C$OMP END MASTER
587C$OMP BARRIER
588
589          lastcff=0.
590          do iter=iif,iif+irhot-1
591          lastcff=lastcff+((iter+1-iif)/rrhot)*
592     &        weight2(iif+irhot-1,iter)
593          enddo
594         
595          lastcff=1./lastcff
596                               
597# ifdef AGRIF_OBC_SOUTH
598      if (SOUTHERN_EDGE) then
599          j=Jstr
600         
601          do i=Istr,Iend
602          DV_south4(i,iif-1) = DV_south3(i,j,iif-1)
603          enddo
604         
605          if (iif+irhot-1<=nfast) then
606           DV_south6=0.
607          do iter=0,iif-1
608          cff1=weight2(iif+irhot-1,iter)
609            do i=Istr,Iend
610             DV_south6(i) = DV_south6(i) +cff1*DV_south4(i,iter)
611            enddo
612          enddo
613         
614          do iter=iif,iif+irhot-2
615          do i=Istr,Iend
616          DV_south6(i)=DV_south6(i)+((iif+irhot-1-iter)/rrhot)*
617     &        weight2(iif+irhot-1,iter)*DV_south4(i,iif-1)
618          enddo
619          enddo
620         
621          do i=Istr,Iend
622          DV_south6(i)=lastcff*((dinterp(i,j)/rrhox)-DV_south6(i))
623#  ifdef MASKING
624     &                               * vmask(i,j)
625#   endif         
626          enddo     
627         
628          do i=Istr,Iend
629          DV_south6(i)=(DV_south6(i)-DV_south4(i,iif-1))/rrhot
630          enddo 
631          endif           
632                                       
633      endif
634# endif 
635
636# ifdef AGRIF_OBC_NORTH
637      if (NORTHERN_EDGE) then 
638          j=Jend+1
639         
640          do i=Istr,Iend
641          DV_north4(i,iif-1) = DV_north3(i,j,iif-1)
642          enddo
643         
644          if (iif+irhot-1<=nfast) then
645           DV_north6=0.
646          do iter=0,iif-1
647          cff1=weight2(iif+irhot-1,iter)
648            do i=Istr,Iend
649             DV_north6(i) = DV_north6(i) +cff1*DV_north4(i,iter)
650            enddo
651          enddo
652         
653          do iter=iif,iif+irhot-2
654          do i=Istr,Iend
655          DV_north6(i)=DV_north6(i)+((iif+irhot-1-iter)/rrhot)*
656     &        weight2(iif+irhot-1,iter)*DV_north4(i,iif-1)
657          enddo
658          enddo
659         
660
661         
662          do i=Istr,Iend
663          DV_north6(i)=lastcff*((dinterp(i,j)/rrhox)-DV_north6(i))
664#  ifdef MASKING
665     &                               * vmask(i,j)
666#   endif             
667          enddo     
668         
669          do i=Istr,Iend
670          DV_north6(i)=(DV_north6(i)-DV_north4(i,iif-1))/rrhot
671          enddo 
672          endif                 
673                                     
674      endif
675# endif
676
677# ifdef AGRIF_OBC_WEST             
678      if (WESTERN_EDGE) then
679          i = Istr-1
680
681          do j=Jstr,Jend
682          DV_west4(j,iif-1) = DV_west3(i,j,iif-1)
683          enddo
684         
685          if (iif+irhot-1<=nfast) then
686           DV_west6=0.
687          do iter=0,iif-1
688          cff1=weight2(iif+irhot-1,iter)
689            do j=Jstr,Jend
690             DV_west6(j) = DV_west6(j) +cff1*DV_west4(j,iter)
691            enddo
692          enddo
693         
694          do iter=iif,iif+irhot-2
695          do j=Jstr,Jend
696          DV_west6(j)=DV_west6(j)+((iif+irhot-1-iter)/rrhot)*
697     &        weight2(iif+irhot-1,iter)*DV_west4(j,iif-1)
698          enddo
699          enddo
700         
701
702         
703          do j=Jstr,Jend
704          DV_west6(j)=lastcff*((dinterp(i,j)/rrhoy)-DV_west6(j))
705#  ifdef MASKING
706     &                               * vmask(i,j)
707#   endif                     
708          enddo     
709         
710          do j=Jstr,Jend
711          DV_west6(j)=(DV_west6(j)-DV_west4(j,iif-1))/rrhot
712          enddo 
713          endif
714                                       
715      endif
716# endif
717
718# ifdef AGRIF_OBC_EAST       
719      if (EASTERN_EDGE) then
720          i=Iend+1         
721         
722          do j=Jstr,Jend
723          DV_east4(j,iif-1) = DV_east3(i,j,iif-1)
724          enddo
725         
726          if (iif+irhot-1<=nfast) then
727           DV_east6=0.
728          do iter=0,iif-1
729          cff1=weight2(iif+irhot-1,iter)
730            do j=Jstr,Jend
731             DV_east6(j) = DV_east6(j) +cff1*DV_east4(j,iter)
732            enddo
733          enddo
734         
735          do iter=iif,iif+irhot-2
736          do j=Jstr,Jend
737          DV_east6(j)=DV_east6(j)+((iif+irhot-1-iter)/rrhot)*
738     &        weight2(iif+irhot-1,iter)*DV_east4(j,iif-1)
739          enddo
740          enddo
741         
742          do j=Jstr,Jend
743          DV_east6(j)=lastcff*((dinterp(i,j)/rrhoy)-DV_east6(j))
744#  ifdef MASKING
745     &                               * vmask(i,j)
746#   endif         
747          enddo     
748         
749          do j=Jstr,Jend
750          DV_east6(j)=(DV_east6(j)-DV_east4(j,iif-1))/rrhot
751          enddo 
752          endif         
753                                 
754      endif
755# endif
756        V2DTimeindex = parentnbstep
757      endif
758         
759
760                     
761# ifdef AGRIF_OBC_SOUTH
762      if (SOUTHERN_EDGE) then
763          j=Jstr
764         
765          do i=Istr,Iend
766           DV_south4(i,iif-1) = DV_south3(i,j,iif-1)
767          enddo
768           
769          do i=Istr,Iend
770            DV_south(i)=DV_south4(i,iif-1)+DV_south6(i)
771          enddo
772                                       
773      endif
774# endif 
775
776# ifdef AGRIF_OBC_NORTH
777      if (NORTHERN_EDGE) then 
778          j=Jend+1
779
780          do i=Istr,Iend
781           DV_north4(i,iif-1) = DV_north3(i,j,iif-1)
782          enddo
783           
784          do i=Istr,Iend
785            DV_north(i)=DV_north4(i,iif-1)+DV_north6(i)
786          enddo
787                                               
788      endif
789# endif
790
791# ifdef AGRIF_OBC_WEST             
792      if (WESTERN_EDGE) then
793          i = Istr-1         
794
795          do j=Jstr,Jend
796           DV_west4(j,iif-1) = DV_west3(i,j,iif-1)
797          enddo
798           
799          do j=Jstr,Jend
800            DV_west(j)=DV_west4(j,iif-1)+DV_west6(j)
801          enddo
802                                       
803      endif
804# endif
805
806# ifdef AGRIF_OBC_EAST       
807      if (EASTERN_EDGE) then
808          i=Iend+1
809         
810          do j=Jstr,Jend
811           DV_east4(j,iif-1) = DV_east3(i,j,iif-1)
812          enddo
813           
814          do j=Jstr,Jend
815            DV_east(j)=DV_east4(j,iif-1)+DV_east6(j)
816          enddo         
817                                 
818      endif
819# endif     
820         
821# if defined M2_FRC_BRY || defined M2NUDGING
822!
823! Apply the value to vbclm or vbarbry
824!
825      cffy = g*dtfast*2./(1.+rrhoy)
826
827#  ifdef AGRIF_2WAY
828      cffy = 0.
829#  endif
830
831#  ifdef AGRIF_OBC_SOUTH
832      if (SOUTHERN_EDGE) then
833         do i=Istr,Iend
834#   ifdef M2_FRC_BRY
835           vbarbry_south(i)=
836#   else
837           vbclm(i,Jstr)=(cffy/on_v(i,Jstr))*
838     &                   (SSH(i,Jstr)-zeta(i,Jstr,knew))   +
839#   endif
840#   ifdef AGRIF_FLUX_BC
841     &              (2.*DV_south(i)/((h(i,Jstr-1)+zeta(i,Jstr-1,knew)
842     &                                  +h(i,Jstr)+zeta(i,Jstr,knew))
843     &                                                 *om_v(i,Jstr)))
844#   else
845     &                  DV_south(i)
846#   endif
847#   ifdef MASKING
848     &                  *vmask(i,Jstr)
849#   endif
850       enddo
851      endif
852#  endif 
853
854#  ifdef AGRIF_OBC_NORTH
855      if (NORTHERN_EDGE) then
856        do i=Istr,Iend
857#   ifdef M2_FRC_BRY
858          vbarbry_north(i)=
859#   else
860          vbclm(i,Jend+1)=-(cffy/on_v(i,Jend+1))*
861     &                     (SSH(i,Jend)-zeta(i,Jend,knew))  +
862#   endif
863#   ifdef AGRIF_FLUX_BC
864     &                  (2.*DV_north(i)/(( h(i,Jend)+zeta(i,Jend,knew)
865     &                               +h(i,Jend+1)+zeta(i,Jend+1,knew))
866     &                                               *om_v(i,Jend+1)))
867#   else
868     &                  DV_north(i)
869#   endif
870#   ifdef MASKING
871     &                  *vmask(i,Jend+1)
872#   endif
873        enddo
874      endif
875#  endif
876
877#  ifdef AGRIF_OBC_EAST       
878      if (EASTERN_EDGE) then
879        do j=Jstr,Jend
880#   ifdef M2_FRC_BRY
881          vbarbry_east(j)=
882#   else
883          vbclm(Iend+1,j)=
884#   endif
885#   ifdef AGRIF_FLUX_BC
886     &             (2.*DV_east(j)/(( h(Iend+1,j)+zeta(Iend+1,j,knew)
887     &                           +h(Iend+1,j-1)+zeta(Iend+1,j-1,knew))
888     &                                             *om_v(Iend+1,j)))
889#   else
890     &                 DV_east(j)
891#   endif
892#   ifdef MASKING
893     &                 *vmask(Iend+1,j)
894#   endif
895        enddo
896      endif       
897#  endif 
898
899#  ifdef AGRIF_OBC_WEST             
900      if (WESTERN_EDGE) then
901        do j=Jstr,Jend
902#   ifdef M2_FRC_BRY
903          vbarbry_west(j)=
904#   else
905          vbclm(Istr-1,j)=
906#   endif
907#   ifdef AGRIF_FLUX_BC
908     &                (2.*DV_west(j)/((h(Istr-1,j)+zeta(Istr-1,j,knew)
909     &                           +h(Istr-1,j-1)+zeta(Istr-1,j-1,knew))
910     &                                                *om_v(Istr-1,j)))
911#   else
912     &                  DV_west(j)
913#   endif
914#   ifdef MASKING
915     &                  *vmask(Istr-1,j)
916#   endif
917        enddo
918      endif
919#  endif
920# endif /* M2_FRC_BRY || M2NUDGING */
921      return
922      end
923
924      subroutine v2Dinterp(tabres,i1,i2,j1,j2)
925      implicit none
926# include "param.h"
927# include "grid.h"
928# include "ocean2d.h"
929# include "scalars.h"
930
931       integer i1,i2,j1,j2
932       real tabres(i1:i2,j1:j2)
933       integer i,j,iter, isize
934       real :: cff1
935!$AGRIF_DO_NOT_TREAT
936      real,dimension(:,:,:),pointer :: coarsevalues
937      integer :: nbgrid, indinterp
938      common/interp2d/coarsevalues, nbgrid, indinterp 
939!$AGRIF_END_DO_NOT_TREAT 
940
941       integer :: oldindinterp
942       
943       isize = (j2-j1+1)*(i2-i1+1)
944       
945       IF (iif == 1) THEN
946       IF (.NOT.Associated(coarsevalues)) THEN
947       Allocate(coarsevalues(isize,0:nfast,1))
948       ELSE       
949       CALL checksizeinterp(indinterp+isize,nfast)
950       ENDIF
951       
952       oldindinterp = indinterp
953       
954       do j=max(j1,lbound(h,2)+1),j2
955       do i=i1,i2
956       oldindinterp=oldindinterp+1
957       coarsevalues(oldindinterp,0,nbgrid) = 
958     &  0.5*(h(i,j-1)+zeta(i,j-1,kstp)+h(i,j)+
959     &  zeta(i,j,kstp))*om_v(i,j)*vbar(i,j,kstp)
960       enddo
961       enddo
962             
963       ENDIF
964       
965       oldindinterp = indinterp
966       
967       do j=max(j1,lbound(h,2)+1),j2
968       do i=i1,i2
969       oldindinterp=oldindinterp+1
970       coarsevalues(oldindinterp,iif,nbgrid) =
971     &   0.5*(h(i,j-1)+zeta(i,j-1,knew)+h(i,j)+
972     &  zeta(i,j,knew))*om_v(i,j)*vbar(i,j,knew)
973       enddo
974       enddo
975         
976       do iter=0,iif
977          cff1 = weight2(iif,iter)
978          call copy1d(tabres,coarsevalues(indinterp+1,iter,nbgrid),
979     &                cff1,isize)
980       enddo
981       
982       indinterp = oldindinterp
983               
984      return
985      end   
986     
987      subroutine copy1d(tab1,tab2,cff,isize)
988      real,dimension(isize) :: tab1,tab2
989      integer :: i
990     
991      do i=1,isize
992      tab1(i)=tab1(i)+cff*tab2(i)
993      enddo
994      return
995      end subroutine copy1d
996
997!
998      subroutine checksizeinterp(isize,nfast)
999      integer :: isize,nfast
1000      real,dimension(:,:,:),allocatable :: tempvalues
1001      integer :: n1,n3
1002!$AGRIF_DO_NOT_TREAT
1003      real,dimension(:,:,:),pointer :: coarsevalues
1004      integer :: nbgrid, indinterp
1005      common/interp2d/coarsevalues, nbgrid, indinterp
1006!$AGRIF_END_DO_NOT_TREAT
1007
1008       IF (size(coarsevalues,1).LT.(isize)) THEN
1009       n1 = size(coarsevalues,1)
1010       n3 = size(coarsevalues,3)
1011       allocate(tempvalues(n1,0:nfast,n3))
1012       tempvalues=coarsevalues(1:n1,0:nfast,1:n3)
1013       deallocate(coarsevalues)
1014       allocate(coarsevalues(isize,0:nfast,n3))
1015       coarsevalues(1:n1,0:nfast,1:n3) = tempvalues
1016             
1017       deallocate(tempvalues)
1018       ELSE IF (nbgrid.GT.size(coarsevalues,3)) THEN
1019       n1 = size(coarsevalues,1)
1020       n3 = size(coarsevalues,3)
1021       allocate(tempvalues(n1,0:nfast,n3))
1022       tempvalues=coarsevalues(1:n1,0:nfast,1:n3)
1023       deallocate(coarsevalues)
1024       allocate(coarsevalues(n1,0:nfast,nbgrid))
1025       coarsevalues(1:n1,0:nfast,1:n3) = tempvalues
1026             
1027       deallocate(tempvalues)       
1028       ENDIF     
1029      return
1030      end                 
1031!
1032!====================================================================
1033!                   subroutine Agrif_zetabc_interp_tile
1034!====================================================================
1035!
1036      subroutine zetabc_interp_tile(Istr,Iend,Jstr,Jend
1037     &        ,Zeta_west4,Zeta_east4,Zeta_west6,Zeta_east6
1038     &        ,Zeta_south4,Zeta_north4,Zeta_south6,Zeta_north6)
1039      use AGRIF_Util
1040!     implicit none
1041# include "param.h"
1042# include "boundary.h"
1043# include "climat.h"
1044# include "grid.h"
1045# include "ocean2d.h"
1046# include "scalars.h"
1047# include "zoom.h"
1048      integer Istr,Iend,Jstr,Jend, i,j, i1, j1
1049      real tinterp
1050      Integer itrcind
1051      INTEGER :: parentknew, parentkstp,nbstep3dparent
1052      real t1,t2,t3,t4,t5,t6,t9,t10,t11,t7
1053      real tin(2), tout(2)
1054      real cff1
1055      real zeta1,zeta2
1056      external zetainterp
1057      integer :: irhot, irhox, irhoy
1058      real :: rrhot
1059      integer :: iter
1060      real :: onemtinterp
1061      integer :: parentnbstep
1062      real,dimension(Istr-1:Istr,PRIVATE_1DETA_SCRATCH_ARRAY) ::
1063     &                    Zeta_west6
1064      real,dimension(Iend:Iend+1,PRIVATE_1DETA_SCRATCH_ARRAY) ::
1065     &                    Zeta_east6     
1066      real,dimension(PRIVATE_1DXI_SCRATCH_ARRAY,Jstr-1:Jstr) ::
1067     &                    Zeta_south6
1068      real,dimension(PRIVATE_1DXI_SCRATCH_ARRAY,Jend:Jend+1) ::
1069     &         Zeta_north6 
1070     
1071      real,dimension(Istr-1:Istr,
1072     &          PRIVATE_1DETA_SCRATCH_ARRAY,0:NWEIGHT) :: Zeta_west4
1073      real,dimension(Iend:Iend+1,
1074     &          PRIVATE_1DETA_SCRATCH_ARRAY,0:NWEIGHT) :: Zeta_east4     
1075      real,dimension(PRIVATE_1DXI_SCRATCH_ARRAY,
1076     &          Jstr-1:Jstr,0:NWEIGHT) :: Zeta_south4
1077      real,dimension(PRIVATE_1DXI_SCRATCH_ARRAY,
1078     &          Jend:Jend+1,0:NWEIGHT) :: Zeta_north4
1079             
1080      real lastcff
1081!$AGRIF_DO_NOT_TREAT
1082      real,dimension(:,:,:),pointer :: coarsevalues
1083      integer :: nbgrid, indinterp
1084      common/interp2d/coarsevalues, nbgrid, indinterp
1085!$AGRIF_END_DO_NOT_TREAT
1086                 
1087#  ifdef MPI
1088#   define LOCALLM Lmmpi
1089#   define LOCALMM Mmmpi
1090#  else
1091#   define LOCALLM Lm
1092#   define LOCALMM Mm
1093#  endif
1094!
1095# include "compute_auxiliary_bounds.h"
1096
1097!
1098!       return     
1099       irhot = Agrif_Irhot()
1100       irhox = Agrif_Irhox()
1101       irhoy = Agrif_Irhoy()
1102       rrhot = real(irhot)
1103
1104
1105
1106      parentnbstep=Agrif_Parent_Nb_Step()
1107
1108      if (ZetaTimeindex .NE. parentnbstep) then
1109     
1110        do J=JstrR,JendR
1111        do i=IstrR,IendR
1112        dinterp(i,j) = 0.
1113        enddo
1114        enddo
1115
1116C$OMP BARRIER
1117C$OMP MASTER
1118
1119            tinterp=1.
1120                         
1121#ifdef MASKING           
1122            Agrif_UseSpecialvalue=.true.
1123#endif           
1124            Agrif_Specialvalue=0.
1125
1126        nbgrid = Agrif_Fixed()
1127        indinterp = 0
1128               
1129            Call Agrif_Bc_variable(dinterp,Zt_avg1,
1130     &                                        calledweight=tinterp,
1131     &                                        procname=zetainterp)
1132            Agrif_UseSpecialvalue=.false. 
1133C$OMP END MASTER
1134C$OMP BARRIER
1135
1136          lastcff=0.
1137          do iter=iif,iif+irhot-1
1138          lastcff=lastcff+((iter+1-iif)/rrhot)*
1139     &        weight2(iif+irhot-1,iter)
1140          enddo
1141         
1142          lastcff=1./lastcff
1143             
1144# ifdef AGRIF_OBC_SOUTH
1145      if (SOUTHERN_EDGE) then
1146
1147          do j=Jstr-1,Jstr
1148          do i=Istr,Iend
1149          Zeta_south4(i,j,iif-1) = Zeta_south3(i,j,iif-1)
1150          enddo
1151          enddo
1152         
1153          if (iif+irhot-1<=nfast) then
1154           Zeta_south6=0.
1155          do iter=0,iif-1
1156          cff1=weight2(iif+irhot-1,iter)
1157            do j=Jstr-1,Jstr
1158            do i=Istr,Iend
1159             Zeta_south6(i,j) = Zeta_south6(i,j)
1160     &                         +cff1*Zeta_south4(i,j,iter)
1161            enddo
1162            enddo
1163          enddo
1164         
1165          do iter=iif,iif+irhot-2
1166          do j=Jstr-1,Jstr         
1167          do i=Istr,Iend
1168          Zeta_south6(i,j)=Zeta_south6(i,j)+((iif+irhot-1-iter)/rrhot)*
1169     &        weight2(iif+irhot-1,iter)*Zeta_south4(i,j,iif-1)
1170          enddo
1171          enddo
1172          enddo         
1173         
1174
1175          do j=Jstr-1,Jstr         
1176          do i=Istr,Iend
1177          Zeta_south6(i,j)=lastcff*(dinterp(i,j)-Zeta_south6(i,j))
1178          enddo
1179          enddo               
1180         
1181          do j=Jstr-1,Jstr
1182          do i=Istr,Iend
1183          Zeta_south6(i,j)=(Zeta_south6(i,j)-Zeta_south4(i,j,iif-1))
1184     &                           /rrhot
1185          enddo 
1186          enddo         
1187          endif
1188                                         
1189      endif
1190# endif 
1191
1192# ifdef AGRIF_OBC_NORTH
1193      if (NORTHERN_EDGE) then
1194             
1195          do j=Jend,Jend+1
1196          do i=Istr,Iend
1197          Zeta_north4(i,j,iif-1) = Zeta_north3(i,j,iif-1)
1198          enddo
1199          enddo
1200         
1201          if (iif+irhot-1<=nfast) then
1202           Zeta_north6=0.
1203          do iter=0,iif-1
1204          cff1=weight2(iif+irhot-1,iter)
1205            do j=Jend,Jend+1
1206            do i=Istr,Iend
1207             Zeta_north6(i,j) = Zeta_north6(i,j)
1208     &                         +cff1*Zeta_north4(i,j,iter)
1209            enddo
1210            enddo
1211          enddo
1212         
1213          do iter=iif,iif+irhot-2
1214          do j=Jend,Jend+1         
1215          do i=Istr,Iend
1216          Zeta_north6(i,j)=Zeta_north6(i,j)+((iif+irhot-1-iter)/rrhot)*
1217     &        weight2(iif+irhot-1,iter)*Zeta_north4(i,j,iif-1)
1218          enddo
1219          enddo
1220          enddo         
1221         
1222
1223          do j=Jend,Jend+1         
1224          do i=Istr,Iend
1225          Zeta_north6(i,j)=lastcff*(dinterp(i,j)-Zeta_north6(i,j))
1226          enddo
1227          enddo               
1228         
1229          do j=Jend,Jend+1
1230          do i=Istr,Iend
1231          Zeta_north6(i,j)=(Zeta_north6(i,j)-Zeta_north4(i,j,iif-1))
1232     &                           /rrhot
1233          enddo 
1234          enddo         
1235          endif                             
1236      endif
1237# endif
1238
1239# ifdef AGRIF_OBC_WEST             
1240      if (WESTERN_EDGE) then           
1241
1242          do j=Jstr,Jend
1243          do i=Istr-1,Istr
1244          Zeta_west4(i,j,iif-1) = Zeta_west3(i,j,iif-1)
1245          enddo
1246          enddo
1247         
1248          if (iif+irhot-1<=nfast) then
1249           Zeta_west6=0.
1250          do iter=0,iif-1
1251          cff1=weight2(iif+irhot-1,iter)
1252            do j=Jstr,Jend
1253            do i=Istr-1,Istr
1254             Zeta_west6(i,j) = Zeta_west6(i,j)
1255     &                         +cff1*Zeta_west4(i,j,iter)
1256            enddo
1257            enddo
1258          enddo
1259         
1260          do iter=iif,iif+irhot-2
1261          do j=Jstr,Jend         
1262          do i=Istr-1,Istr
1263          Zeta_west6(i,j)=Zeta_west6(i,j)+((iif+irhot-1-iter)/rrhot)*
1264     &        weight2(iif+irhot-1,iter)*Zeta_west4(i,j,iif-1)
1265          enddo
1266          enddo
1267          enddo         
1268         
1269
1270          do j=Jstr,Jend         
1271          do i=Istr-1,Istr
1272          Zeta_west6(i,j)=lastcff*(dinterp(i,j)-Zeta_west6(i,j))
1273          enddo
1274          enddo               
1275         
1276          do j=Jstr,Jend
1277          do i=Istr-1,Istr
1278          Zeta_west6(i,j)=(Zeta_west6(i,j)-Zeta_west4(i,j,iif-1))
1279     &                           /rrhot
1280          enddo 
1281          enddo         
1282          endif
1283                                         
1284      endif
1285# endif
1286
1287# ifdef AGRIF_OBC_EAST       
1288      if (EASTERN_EDGE) then
1289
1290          do j=Jstr,Jend
1291          do i=Iend,Iend+1
1292          Zeta_east4(i,j,iif-1) = Zeta_east3(i,j,iif-1)
1293          enddo
1294          enddo
1295         
1296          if (iif+irhot-1<=nfast) then
1297           Zeta_east6=0.
1298          do iter=0,iif-1
1299          cff1=weight2(iif+irhot-1,iter)
1300            do j=Jstr,Jend
1301            do i=Iend,Iend+1
1302             Zeta_east6(i,j) = Zeta_east6(i,j)
1303     &                         +cff1*Zeta_east4(i,j,iter)
1304            enddo
1305            enddo
1306          enddo
1307         
1308          do iter=iif,iif+irhot-2
1309          do j=Jstr,Jend         
1310          do i=Iend,Iend+1
1311          Zeta_east6(i,j)=Zeta_east6(i,j)+((iif+irhot-1-iter)/rrhot)*
1312     &        weight2(iif+irhot-1,iter)*Zeta_east4(i,j,iif-1)
1313          enddo
1314          enddo
1315          enddo         
1316         
1317
1318          do j=Jstr,Jend         
1319          do i=Iend,Iend+1
1320          Zeta_east6(i,j)=lastcff*(dinterp(i,j)-Zeta_east6(i,j))
1321          enddo
1322          enddo               
1323         
1324          do j=Jstr,Jend
1325          do i=Iend,Iend+1
1326          Zeta_east6(i,j)=(Zeta_east6(i,j)-Zeta_east4(i,j,iif-1))
1327     &                           /rrhot
1328          enddo 
1329          enddo         
1330          endif
1331                                         
1332      endif
1333# endif
1334            ZetaTimeindex = parentnbstep
1335          endif
1336
1337#  ifdef AGRIF_OBC_SOUTH             
1338          if (SOUTHERN_EDGE) then
1339         
1340          do j=Jstr-1,Jstr
1341          do i=Istr,Iend
1342          Zeta_south4(i,j,iif-1) = Zeta_south3(i,j,iif-1)
1343          enddo
1344          enddo
1345#   ifdef ZNUDGING
1346          do j=Jstr-1,Jstr
1347          do i=Istr,Iend
1348          SSH(i,j)=Zeta_south4(i,j,iif-1)+Zeta_south6(i,j)
1349          enddo
1350          enddo
1351#   endif
1352         
1353          endif
1354#  endif
1355#  ifdef AGRIF_OBC_NORTH 
1356          if (NORTHERN_EDGE) then
1357         
1358          do j=Jend,Jend+1
1359          do i=Istr,Iend
1360          Zeta_north4(i,j,iif-1) = Zeta_north3(i,j,iif-1)
1361          enddo
1362          enddo
1363#   ifdef ZNUDGING
1364          do j=Jend,Jend+1
1365          do i=Istr,Iend
1366          SSH(i,j)=Zeta_north4(i,j,iif-1)+Zeta_north6(i,j)
1367          enddo
1368          enddo
1369#   endif
1370          endif
1371#  endif
1372#  ifdef AGRIF_OBC_WEST 
1373          if (WESTERN_EDGE) then
1374         
1375          do j=Jstr,Jend
1376          do i=Istr-1,Istr
1377          Zeta_west4(i,j,iif-1) = Zeta_west3(i,j,iif-1)
1378          enddo
1379          enddo
1380#   ifdef ZNUDGING
1381          do j=Jstr,Jend
1382          do i=Istr-1,Istr
1383          SSH(i,j)=Zeta_west4(i,j,iif-1)+Zeta_west6(i,j)
1384          enddo
1385          enddo
1386#   endif
1387          endif
1388#  endif
1389#  ifdef AGRIF_OBC_EAST
1390          if (EASTERN_EDGE) then
1391         
1392          do j=Jstr,Jend
1393          do i=Iend,Iend+1
1394          Zeta_east4(i,j,iif-1) = Zeta_east3(i,j,iif-1)
1395          enddo
1396          enddo
1397#   ifdef ZNUDGING
1398          do j=Jstr,Jend
1399          do i=Iend,Iend+1
1400          SSH(i,j)=Zeta_east4(i,j,iif-1)+Zeta_east6(i,j)
1401          enddo
1402          enddo
1403#   endif
1404          endif
1405#  endif
1406     
1407      return
1408      end
1409     
1410      subroutine zetainterp(tabres,i1,i2,j1,j2)
1411      implicit none
1412# include "param.h"
1413# include "ocean2d.h"
1414# include "scalars.h"
1415
1416       integer i1,i2,j1,j2
1417       real tabres(i1:i2,j1:j2)
1418       integer i,j,iter, isize
1419       real :: cff1
1420!$AGRIF_DO_NOT_TREAT
1421      real,dimension(:,:,:),pointer :: coarsevalues
1422      integer :: nbgrid, indinterp
1423      common/interp2d/coarsevalues, nbgrid, indinterp 
1424!$AGRIF_END_DO_NOT_TREAT 
1425
1426       integer :: oldindinterp
1427       
1428       isize = (j2-j1+1)*(i2-i1+1)
1429       
1430       IF (iif == 1) THEN
1431       IF (.NOT.Associated(coarsevalues)) THEN
1432       Allocate(coarsevalues(isize,0:nfast,1))
1433       ELSE       
1434       CALL checksizeinterp(indinterp+isize,nfast)
1435       ENDIF
1436       
1437       oldindinterp = indinterp
1438       
1439       do j=j1,j2
1440       do i=i1,i2
1441       oldindinterp=oldindinterp+1
1442       coarsevalues(oldindinterp,0,nbgrid) = 
1443     &  zeta(i,j,kstp)
1444       enddo
1445       enddo
1446             
1447       ENDIF
1448       
1449       oldindinterp = indinterp
1450       
1451       do j=j1,j2
1452       do i=i1,i2
1453       oldindinterp=oldindinterp+1
1454       coarsevalues(oldindinterp,iif,nbgrid) =
1455     &   zeta(i,j,knew)
1456       enddo
1457       enddo
1458         
1459       do iter=0,iif
1460          cff1 = weight2(iif,iter)
1461          call copy1d(tabres,coarsevalues(indinterp+1,iter,nbgrid),
1462     &                cff1,isize)
1463       enddo
1464       
1465       indinterp = oldindinterp
1466     
1467      return
1468      end     
1469#else
1470      subroutine zommbc_2D_empty()
1471      return
1472      end
1473#endif
Note: See TracBrowser for help on using the repository browser.