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_user.F90 in utils/tools/DOMAINcfg/src – NEMO

source: utils/tools/DOMAINcfg/src/agrif_user.F90 @ 14698

Last change on this file since 14698 was 14698, checked in by jchanut, 3 years ago

#2638: initialize agrif ghost cells number (to zero) on root grid (DOMAINcfg)

File size: 45.9 KB
Line 
1#if defined key_agrif
2   SUBROUTINE agrif_user()
3      !!----------------------------------------------------------------------
4      !!                 *** ROUTINE agrif_user ***
5      !!----------------------------------------------------------------------
6   END SUBROUTINE agrif_user
7
8   SUBROUTINE agrif_initworkspace()
9      !!----------------------------------------------------------------------
10      !!                 *** ROUTINE Agrif_InitWorkspace ***
11      !!----------------------------------------------------------------------
12   END SUBROUTINE agrif_initworkspace
13
14   SUBROUTINE Agrif_InitValues
15      !!----------------------------------------------------------------------
16      !!                 *** ROUTINE Agrif_InitValues ***
17      !!
18      !! ** Purpose :: Declaration of variables to be interpolated
19      !!----------------------------------------------------------------------
20      USE Agrif_Util
21      USE dom_oce
22      USE nemogcm
23      USE domain
24      !!
25      IMPLICIT NONE
26
27      ! No temporal refinement
28      CALL Agrif_Set_coeffreft(1)
29
30      CALL nemo_init       !* Initializations of each fine grid
31
32      CALL dom_nam
33
34   END SUBROUTINE Agrif_InitValues
35
36   SUBROUTINE Agrif_InitValues_cont
37      !!----------------------------------------------------------------------
38      !!                 *** ROUTINE Agrif_InitValues_cont ***
39      !!
40      !! ** Purpose :: Initialisation of variables to be interpolated
41      !!----------------------------------------------------------------------
42      USE dom_oce
43      USE lbclnk
44      !!
45      IMPLICIT NONE
46      !
47      INTEGER :: irafx, irafy
48      LOGICAL :: ln_perio, l_deg
49      !
50      irafx = agrif_irhox()
51      irafy = agrif_irhoy()
52
53
54   !       IF(jperio /=1 .AND. jperio/=4 .AND. jperio/=6 ) THEN
55   !          nx = (nbcellsx)+2*nbghostcellsfine+2
56   !          ny = (nbcellsy)+2*nbghostcellsfine+2
57   !          nbghostcellsfine_tot_x= nbghostcellsfine_x +1
58   !          nbghostcellsfine_tot_y= nbghostcellsfine_y +1
59   !       ELSE
60   !         nx = (nbcellsx)+2*nbghostcellsfine_x
61   !         ny = (nbcellsy)+2*nbghostcellsfine+2
62   !         nbghostcellsfine_tot_x= 1
63   !         nbghostcellsfine_tot_y= nbghostcellsfine_y +1
64   !      ENDIF
65   !    ELSE
66   !       nbghostcellsfine = 0
67   !       nx = nbcellsx+irafx
68   !       ny = nbcellsy+irafy
69
70      WRITE(*,*) ' '
71      WRITE(*,*)'Size of the High resolution grid: ',jpi,' x ',jpj
72      WRITE(*,*) ' '
73      ln_perio = .FALSE.
74      l_deg = .TRUE.
75 
76      IF( jperio == 1 .OR. jperio == 2 .OR. jperio == 4 ) ln_perio=.TRUE.
77      IF ( Agrif_Parent(jphgr_msh)==2 &
78      &.OR.Agrif_Parent(jphgr_msh)==3 & 
79      &.OR.Agrif_Parent(jphgr_msh)==5 ) l_deg = .FALSE.
80
81      CALL agrif_init_lonlat()
82   
83      IF ( l_deg ) THEN
84         WHERE (glamt < -180) glamt = glamt +360.
85         WHERE (glamt > +180) glamt = glamt -360.
86         WHERE (glamu < -180) glamu = glamu +360.
87         WHERE (glamu > +180) glamu = glamu -360.
88         WHERE (glamv < -180) glamv = glamv +360.
89         WHERE (glamv > +180) glamv = glamv -360.
90         WHERE (glamf < -180) glamf = glamf +360.
91         WHERE (glamf > +180) glamf = glamf -360.
92      ENDIF
93
94      CALL lbc_lnk( 'glamt', glamt, 'T', 1._wp)
95      CALL lbc_lnk( 'gphit', gphit, 'T', 1._wp)
96
97      CALL lbc_lnk( 'glamu', glamu, 'U', 1._wp)
98      CALL lbc_lnk( 'gphiu', gphiu, 'U', 1._wp)
99
100      CALL lbc_lnk( 'glamv', glamv, 'V', 1._wp)
101      CALL lbc_lnk( 'gphiv', gphiv, 'V', 1._wp)
102
103      CALL lbc_lnk( 'glamf', glamf, 'F', 1._wp)
104      CALL lbc_lnk( 'gphif', gphif, 'F', 1._wp)
105
106      ! Correct South and North
107      IF ((nbondj == -1).OR.(nbondj == 2)) THEN
108         glamt(:,1) = glamt(:,2)
109         gphit(:,1) = gphit(:,2)
110         glamu(:,1) = glamu(:,2)
111         gphiu(:,1) = gphiu(:,2)
112         glamv(:,1) = glamv(:,2)
113         gphiv(:,1) = gphiv(:,2)
114      ENDIF
115      IF ((nbondj == 1).OR.(nbondj == 2)) THEN
116         glamt(:,jpj) = glamt(:,jpj-1)
117         gphit(:,jpj) = gphit(:,jpj-1)
118         glamu(:,jpj) = glamu(:,jpj-1)
119         gphiu(:,jpj) = gphiu(:,jpj-1)
120         glamv(:,jpj) = glamv(:,jpj-1)
121         gphiv(:,jpj) = gphiv(:,jpj-1)
122         glamf(:,jpj) = glamf(:,jpj-1)
123         gphif(:,jpj) = gphif(:,jpj-1)
124      ENDIF
125
126      ! Correct West and East
127      IF( jperio /= 1 ) THEN
128         IF((nbondi == -1) .OR. (nbondi == 2) ) THEN
129            glamt(1,:) = glamt(2,:)
130            gphit(1,:) = gphit(2,:)
131            glamu(1,:) = glamu(2,:)
132            gphiu(1,:) = gphiu(2,:)
133            glamv(1,:) = glamv(2,:)
134            gphiv(1,:) = gphiv(2,:)
135         ENDIF
136         IF( (nbondi == 1) .OR. (nbondi == 2) ) THEN
137            glamt(jpi,:) = glamt(jpi-1,:)
138            gphit(jpi,:) = gphit(jpi-1,:)
139            glamu(jpi,:) = glamu(jpi-1,:)
140            gphiu(jpi,:) = gphiu(jpi-1,:)
141            glamv(jpi,:) = glamv(jpi-1,:)
142            gphiv(jpi,:) = gphiv(jpi-1,:)
143            glamf(jpi,:) = glamf(jpi-1,:)
144            gphif(jpi,:) = gphif(jpi-1,:)
145         ENDIF
146      ENDIF
147
148      CALL agrif_init_scales()
149
150      ! Correct South and North
151      IF( (nbondj == -1) .OR. (nbondj == 2) ) THEN
152         e1t(:,1) = e1t(:,2)
153         e2t(:,1) = e2t(:,2)
154         e1u(:,1) = e1u(:,2)
155         e2u(:,1) = e2u(:,2)
156         e1v(:,1) = e1v(:,2)
157         e2v(:,1) = e2v(:,2)
158      ENDIF
159      IF( (nbondj == 1) .OR. (nbondj == 2) ) THEN
160         e1t(:,jpj) = e1t(:,jpj-1)
161         e2t(:,jpj) = e2t(:,jpj-1)
162         e1u(:,jpj) = e1u(:,jpj-1)
163         e2u(:,jpj) = e2u(:,jpj-1)
164         e1v(:,jpj) = e1v(:,jpj-1)
165         e2v(:,jpj) = e2v(:,jpj-1)
166         e1f(:,jpj) = e1f(:,jpj-1)
167         e2f(:,jpj) = e2f(:,jpj-1)
168      ENDIF
169
170      ! Correct West and East
171      IF( jperio /= 1 ) THEN
172         IF( (nbondj == -1) .OR. (nbondj == 2) ) THEN
173            e1t(1,:) = e1t(2,:)
174            e2t(1,:) = e2t(2,:)
175            e1u(1,:) = e1u(2,:)
176            e2u(1,:) = e2u(2,:)
177            e1v(1,:) = e1v(2,:)
178            e2v(1,:) = e2v(2,:)
179         ENDIF
180         IF( (nbondj == 1) .OR. (nbondj == 2) ) THEN
181            e1t(jpi,:) = e1t(jpi-1,:)
182            e2t(jpi,:) = e2t(jpi-1,:)
183            e1u(jpi,:) = e1u(jpi-1,:)
184            e2u(jpi,:) = e2u(jpi-1,:)
185            e1v(jpi,:) = e1v(jpi-1,:)
186            e2v(jpi,:) = e2v(jpi-1,:)
187            e1f(jpi,:) = e1f(jpi-1,:)
188            e2f(jpi,:) = e2f(jpi-1,:)
189         ENDIF
190      ENDIF
191
192   END SUBROUTINE Agrif_InitValues_cont
193
194
195   SUBROUTINE agrif_declare_var()
196      !!----------------------------------------------------------------------
197      !!                 *** ROUTINE Agrif_InitValues_cont ***
198      !!
199      !! ** Purpose :: Declaration of variables to be interpolated
200      !!----------------------------------------------------------------------
201      USE par_oce
202      USE dom_oce
203      USE agrif_profiles
204      USE agrif_parameters
205
206      IMPLICIT NONE
207
208      INTEGER :: ind1, ind2, ind3
209      INTEGER ::nbghostcellsfine_tot_x, nbghostcellsfine_tot_y
210      INTEGER :: irafx
211
212      EXTERNAL :: nemo_mapping
213
214      ! 1. Declaration of the type of variable which have to be interpolated
215      !---------------------------------------------------------------------
216
217      ind2 = nn_hls + 1 + nbghostcells_x
218      ind3 = nn_hls + 1 + nbghostcells_y_s
219
220      nbghostcellsfine_tot_x=nbghostcells_x+1
221      nbghostcellsfine_tot_y=max(nbghostcells_y_s,nbghostcells_y_n)+1
222
223      irafx = Agrif_irhox()
224
225      ! In case of East-West periodicity, prevent AGRIF interpolation at east and west boundaries
226      ! The procnames will not be CALLed at these boundaries
227      if (jperio == 1) THEN
228        CALL Agrif_Set_NearCommonBorderX(.TRUE.)
229        CALL Agrif_Set_DistantCommonBorderX(.TRUE.)
230      endif
231      if (.not.lk_south) THEN
232        CALL Agrif_Set_NearCommonBorderY(.TRUE.)
233      endif
234
235      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/jpi,jpj/),glamt_id)
236      CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/jpi,jpj/),glamu_id)
237      CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/jpi,jpj/),glamv_id)
238      CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/jpi,jpj/),glamf_id)
239
240      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/jpi,jpj/),gphit_id)
241      CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/jpi,jpj/),gphiu_id)
242      CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/jpi,jpj/),gphiv_id)
243      CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/jpi,jpj/),gphif_id)
244
245      ! Horizontal scale factors
246
247      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/jpi,jpj/),e1t_id)
248      CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/jpi,jpj/),e1u_id)
249      CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/jpi,jpj/),e1v_id)
250      CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/jpi,jpj/),e1f_id)
251
252      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/jpi,jpj/),e2t_id)
253      CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/jpi,jpj/),e2u_id)
254      CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/jpi,jpj/),e2v_id)
255      CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/jpi,jpj/),e2f_id)
256
257      ! Bathymetry
258
259      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/jpi,jpj/),bathy_id)
260
261      ! Vertical scale factors
262      CALL agrif_declare_variable((/2,2,0/),(/ind2,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/jpi,jpj,jpk/),e3t_id)
263      CALL agrif_declare_variable((/2,2,0/),(/ind2,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/jpi,jpj,jpk/),e3t_copy_id)
264      CALL agrif_declare_variable((/2,2,0/),(/ind2,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/jpi,jpj,jpk+1/),e3t_connect_id)
265
266      CALL agrif_declare_variable((/1,2,0/),(/ind2-1,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/jpi,jpj,jpk/),e3u_id)
267      CALL agrif_declare_variable((/2,1,0/),(/ind2,ind3-1,0/),(/'x','y','N'/),(/1,1,1/),(/jpi,jpj,jpk/),e3v_id)
268
269      ! Bottom level
270
271      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/jpi,jpj/),bottom_level_id)
272
273      CALL Agrif_Set_bcinterp(glamt_id,interp=AGRIF_linear)
274      CALL Agrif_Set_interp(glamt_id,interp=AGRIF_linear)
275      CALL Agrif_Set_bc( glamt_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) )
276
277      CALL Agrif_Set_bcinterp(glamu_id,interp=AGRIF_linear)
278      CALL Agrif_Set_interp(glamu_id,interp=AGRIF_linear)
279      CALL Agrif_Set_bc( glamu_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
280
281      CALL Agrif_Set_bcinterp(glamv_id,interp=AGRIF_linear)
282      CALL Agrif_Set_interp(glamv_id,interp=AGRIF_linear)
283      CALL Agrif_Set_bc( glamv_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
284
285      CALL Agrif_Set_bcinterp(glamf_id,interp=AGRIF_linear)
286      CALL Agrif_Set_interp(glamf_id,interp=AGRIF_linear)
287      CALL Agrif_Set_bc( glamf_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
288
289      CALL Agrif_Set_bcinterp(gphit_id,interp=AGRIF_linear)
290      CALL Agrif_Set_interp(gphit_id,interp=AGRIF_linear)
291      CALL Agrif_Set_bc( gphit_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) )
292
293      CALL Agrif_Set_bcinterp(gphiu_id,interp=AGRIF_linear)
294      CALL Agrif_Set_interp(gphiu_id,interp=AGRIF_linear)
295      CALL Agrif_Set_bc( gphiu_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
296
297      CALL Agrif_Set_bcinterp(gphiv_id,interp=AGRIF_linear)
298      CALL Agrif_Set_interp(gphiv_id,interp=AGRIF_linear)
299      CALL Agrif_Set_bc( gphiv_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
300
301      CALL Agrif_Set_bcinterp(gphif_id,interp=AGRIF_linear)
302      CALL Agrif_Set_interp(gphif_id,interp=AGRIF_linear)
303      CALL Agrif_Set_bc( gphif_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
304
305      !
306
307      CALL Agrif_Set_bcinterp(e1t_id,interp=AGRIF_ppm)
308      CALL Agrif_Set_interp(e1t_id,interp=AGRIF_ppm)
309      CALL Agrif_Set_bc( e1t_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) )
310
311      CALL Agrif_Set_bcinterp(e1u_id, interp1=Agrif_linear, interp2=AGRIF_ppm)
312      CALL Agrif_Set_interp(e1u_id, interp1=Agrif_linear, interp2=AGRIF_ppm)
313      CALL Agrif_Set_bc( e1u_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
314
315      CALL Agrif_Set_bcinterp(e1v_id,interp1=AGRIF_ppm, interp2=Agrif_linear)
316      CALL Agrif_Set_interp(e1v_id, interp1=AGRIF_ppm, interp2=Agrif_linear)
317      CALL Agrif_Set_bc( e1v_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
318
319      CALL Agrif_Set_bcinterp(e1f_id,interp=AGRIF_linear)
320      CALL Agrif_Set_interp(e1f_id,interp=AGRIF_linear)
321      CALL Agrif_Set_bc( e1f_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
322
323      CALL Agrif_Set_bcinterp(e2t_id,interp=AGRIF_ppm)
324      CALL Agrif_Set_interp(e2t_id,interp=AGRIF_ppm)
325      CALL Agrif_Set_bc( e2t_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) )
326
327      CALL Agrif_Set_bcinterp(e2u_id,interp1=Agrif_linear, interp2=AGRIF_ppm)
328      CALL Agrif_Set_interp(e2u_id,interp1=Agrif_linear, interp2=AGRIF_ppm)
329      CALL Agrif_Set_bc( e2u_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
330
331      CALL Agrif_Set_bcinterp(e2v_id,interp1=AGRIF_ppm, interp2=Agrif_linear)
332      CALL Agrif_Set_interp(e2v_id,interp1=AGRIF_ppm, interp2=Agrif_linear)
333      CALL Agrif_Set_bc( e2v_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
334
335      CALL Agrif_Set_bcinterp(e2f_id,interp=AGRIF_linear)
336      CALL Agrif_Set_interp(e2f_id,interp=AGRIF_linear)
337      CALL Agrif_Set_bc( e2f_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) )
338
339      CALL Agrif_Set_bcinterp(bathy_id,interp=AGRIF_linear)
340      CALL Agrif_Set_interp(bathy_id,interp=AGRIF_linear)
341      CALL Agrif_Set_bc( bathy_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) )
342
343      ! Vertical scale factors
344      CALL Agrif_Set_bcinterp(e3t_id,interp=AGRIF_ppm)
345      CALL Agrif_Set_interp(e3t_id,interp=AGRIF_ppm)
346      CALL Agrif_Set_bc( e3t_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) )
347      CALL Agrif_Set_Updatetype( e3t_id, update = AGRIF_Update_Average)
348
349      CALL Agrif_Set_bcinterp(e3t_copy_id,interp=AGRIF_constant)
350      CALL Agrif_Set_interp(e3t_copy_id,interp=AGRIF_constant)
351      CALL Agrif_Set_bc( e3t_copy_id, (/-npt_copy*irafx-1,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/))
352
353!      CALL Agrif_Set_bcinterp(e3t_connect_id,interp=AGRIF_linear)
354!      CALL Agrif_Set_interp(e3t_connect_id,interp=AGRIF_linear)
355      CALL Agrif_Set_bcinterp(e3t_connect_id,interp=AGRIF_constant)
356      CALL Agrif_Set_interp(e3t_connect_id,interp=AGRIF_constant)
357      CALL Agrif_Set_bc( e3t_connect_id, (/-(npt_copy+npt_connect)*irafx-1,-npt_copy*irafx-1/))
358
359      CALL Agrif_Set_bcinterp(e3u_id, interp1=Agrif_linear, interp2=AGRIF_ppm)
360      CALL Agrif_Set_interp(e3u_id, interp1=Agrif_linear, interp2=AGRIF_ppm)
361      CALL Agrif_Set_bc( e3u_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) )
362      CALL Agrif_Set_Updatetype(e3u_id,update1 = Agrif_Update_Copy, update2 = Agrif_Update_Average)
363
364      CALL Agrif_Set_bcinterp(e3v_id,interp1=AGRIF_ppm, interp2=Agrif_linear)
365      CALL Agrif_Set_interp(e3v_id, interp1=AGRIF_ppm, interp2=Agrif_linear)
366      CALL Agrif_Set_bc( e3v_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) )
367      CALL Agrif_Set_Updatetype(e3v_id,update1 = Agrif_Update_Average, update2 = Agrif_Update_Copy)
368
369      ! Bottom level
370      CALL Agrif_Set_bcinterp(bottom_level_id,interp=AGRIF_constant)
371      CALL Agrif_Set_interp(bottom_level_id,interp=AGRIF_constant)
372      CALL Agrif_Set_bc( bottom_level_id, (/-npt_copy*irafx-1,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/))
373      CALL Agrif_Set_Updatetype( bottom_level_id, update = AGRIF_Update_Max)
374
375      CALL Agrif_Set_ExternalMapping(nemo_mapping)
376
377   END SUBROUTINE agrif_declare_var
378
379   SUBROUTINE nemo_mapping(ndim,ptx,pty,bounds,bounds_chunks,correction_required,nb_chunks)
380      USE dom_oce
381      INTEGER :: ndim
382      INTEGER :: ptx, pty
383      INTEGER,DIMENSION(ndim,2,2) :: bounds
384      INTEGER,DIMENSION(:,:,:,:),allocatable :: bounds_chunks
385      LOGICAL,DIMENSION(:),allocatable :: correction_required
386      INTEGER :: nb_chunks
387      INTEGER :: i
388
389      IF (agrif_debug_interp) THEN
390         DO i = 1, ndim
391             print *,'direction = ',i,bounds(i,1,2),bounds(i,2,2)
392         END DO
393      ENDIF
394
395      IF( bounds(2,2,2) > jpjglo ) THEN
396         IF( bounds(2,1,2) <= jpjglo ) THEN
397            nb_chunks = 2
398            ALLOCATE(bounds_chunks(nb_chunks,ndim,2,2))
399            ALLOCATE(correction_required(nb_chunks))
400            DO i = 1, nb_chunks
401               bounds_chunks(i,:,:,:) = bounds
402            END DO
403
404         ! FIRST CHUNCK (for j<=jpjglo)
405
406            ! Original indices
407            bounds_chunks(1,1,1,1) = bounds(1,1,2)
408            bounds_chunks(1,1,2,1) = bounds(1,2,2)
409            bounds_chunks(1,2,1,1) = bounds(2,1,2)
410            bounds_chunks(1,2,2,1) = jpjglo
411
412            bounds_chunks(1,1,1,2) = bounds(1,1,2)
413            bounds_chunks(1,1,2,2) = bounds(1,2,2)
414            bounds_chunks(1,2,1,2) = bounds(2,1,2)
415            bounds_chunks(1,2,2,2) = jpjglo
416
417            ! Correction required or not
418            correction_required(1)=.FALSE.
419
420         ! SECOND CHUNCK (for j>jpjglo)
421
422            !Original indices
423            bounds_chunks(2,1,1,1) = bounds(1,1,2)
424            bounds_chunks(2,1,2,1) = bounds(1,2,2)
425            bounds_chunks(2,2,1,1) = jpjglo-2
426            bounds_chunks(2,2,2,1) = bounds(2,2,2)
427
428           ! Where to find them
429           ! We use the relation TAB(ji,jj)=TAB(jpiglo-ji+2,jpjglo-2-(jj-jpjglo))
430
431            IF (ptx == 2) THEN ! T, V points
432               bounds_chunks(2,1,1,2) = jpiglo-bounds(1,2,2)+2
433               bounds_chunks(2,1,2,2) = jpiglo-bounds(1,1,2)+2
434            ELSE ! U, F points
435               bounds_chunks(2,1,1,2) = jpiglo-bounds(1,2,2)+1
436               bounds_chunks(2,1,2,2) = jpiglo-bounds(1,1,2)+1
437            ENDIF
438
439            IF (pty == 2) THEN ! T, U points
440               bounds_chunks(2,2,1,2) = jpjglo-2-(bounds(2,2,2) -jpjglo)
441               bounds_chunks(2,2,2,2) = jpjglo-2-(jpjglo-2      -jpjglo)
442            ELSE ! V, F points
443               bounds_chunks(2,2,1,2) = jpjglo-3-(bounds(2,2,2) -jpjglo)
444               bounds_chunks(2,2,2,2) = jpjglo-3-(jpjglo-2      -jpjglo)
445            ENDIF
446     
447            ! Correction required or not
448            correction_required(2)=.TRUE.
449
450         ELSE
451           
452            nb_chunks = 1
453            ALLOCATE(bounds_chunks(nb_chunks,ndim,2,2))
454            ALLOCATE(correction_required(nb_chunks))
455            DO i=1,nb_chunks
456                bounds_chunks(i,:,:,:) = bounds
457            END DO
458
459            bounds_chunks(1,1,1,1) = bounds(1,1,2)
460            bounds_chunks(1,1,2,1) = bounds(1,2,2)
461            bounds_chunks(1,2,1,1) = bounds(2,1,2)
462            bounds_chunks(1,2,2,1) = bounds(2,2,2)
463
464            bounds_chunks(1,1,1,2) = jpiglo-bounds(1,2,2)+2
465            bounds_chunks(1,1,2,2) = jpiglo-bounds(1,1,2)+2
466
467            bounds_chunks(1,2,1,2) = jpjglo-2-(bounds(2,2,2)-jpjglo)
468            bounds_chunks(1,2,2,2) = jpjglo-2-(bounds(2,1,2)-jpjglo)
469
470            IF (ptx == 2) THEN ! T, V points
471               bounds_chunks(1,1,1,2) = jpiglo-bounds(1,2,2)+2
472               bounds_chunks(1,1,2,2) = jpiglo-bounds(1,1,2)+2
473            ELSE ! U, F points
474               bounds_chunks(1,1,1,2) = jpiglo-bounds(1,2,2)+1
475               bounds_chunks(1,1,2,2) = jpiglo-bounds(1,1,2)+1
476            ENDIF
477
478            IF (pty == 2) THEN ! T, U points
479               bounds_chunks(1,2,1,2) = jpjglo-2-(bounds(2,2,2) -jpjglo)
480               bounds_chunks(1,2,2,2) = jpjglo-2-(bounds(2,1,2) -jpjglo)
481            ELSE ! V, F points
482               bounds_chunks(1,2,1,2) = jpjglo-3-(bounds(2,2,2) -jpjglo)
483               bounds_chunks(1,2,2,2) = jpjglo-3-(bounds(2,1,2) -jpjglo)
484            ENDIF
485
486            correction_required(1)=.TRUE.
487
488         ENDIF  ! bounds(2,1,2) <= jpjglo
489
490      ELSE IF (bounds(1,1,2) < 1) THEN
491         
492         IF (bounds(1,2,2) > 0) THEN
493            nb_chunks = 2
494            ALLOCATE(correction_required(nb_chunks))
495            correction_required=.FALSE.
496            ALLOCATE(bounds_chunks(nb_chunks,ndim,2,2))
497            DO i=1,nb_chunks
498               bounds_chunks(i,:,:,:) = bounds
499            END DO
500
501            bounds_chunks(1,1,1,2) = bounds(1,1,2)+jpiglo-2
502            bounds_chunks(1,1,2,2) = 1+jpiglo-2
503
504            bounds_chunks(1,1,1,1) = bounds(1,1,2)
505            bounds_chunks(1,1,2,1) = 1
506
507            bounds_chunks(2,1,1,2) = 2
508            bounds_chunks(2,1,2,2) = bounds(1,2,2)
509
510            bounds_chunks(2,1,1,1) = 2
511            bounds_chunks(2,1,2,1) = bounds(1,2,2)
512         ELSE
513            nb_chunks = 1
514            ALLOCATE(correction_required(nb_chunks))
515            correction_required=.FALSE.
516            ALLOCATE(bounds_chunks(nb_chunks,ndim,2,2))
517            DO i=1,nb_chunks
518               bounds_chunks(i,:,:,:) = bounds
519            END DO
520            bounds_chunks(1,1,1,2) = bounds(1,1,2)+jpiglo-2
521            bounds_chunks(1,1,2,2) = bounds(1,2,2)+jpiglo-2
522
523            bounds_chunks(1,1,1,1) = bounds(1,1,2)
524            bounds_chunks(1,1,2,1) = bounds(1,2,2)
525         ENDIF
526     
527      ELSE
528     
529         nb_chunks=1
530         ALLOCATE(correction_required(nb_chunks))
531         correction_required=.FALSE.
532         ALLOCATE(bounds_chunks(nb_chunks,ndim,2,2))
533         DO i=1,nb_chunks
534            bounds_chunks(i,:,:,:) = bounds
535         END DO
536         bounds_chunks(1,1,1,2) = bounds(1,1,2)
537         bounds_chunks(1,1,2,2) = bounds(1,2,2)
538         bounds_chunks(1,2,1,2) = bounds(2,1,2)
539         bounds_chunks(1,2,2,2) = bounds(2,2,2)
540
541         bounds_chunks(1,1,1,1) = bounds(1,1,2)
542         bounds_chunks(1,1,2,1) = bounds(1,2,2)
543         bounds_chunks(1,2,1,1) = bounds(2,1,2)
544         bounds_chunks(1,2,2,1) = bounds(2,2,2)
545
546      ENDIF
547
548   END SUBROUTINE nemo_mapping
549
550   FUNCTION agrif_external_switch_index(ptx,pty,i1,isens)
551      USE dom_oce
552      INTEGER :: ptx, pty, i1, isens
553      INTEGER :: agrif_external_switch_index
554
555      IF( isens == 1 )  THEN
556         IF( ptx == 2 ) THEN ! T, V points
557            agrif_external_switch_index = jpiglo-i1+2
558         ELSE ! U, F points
559            agrif_external_switch_index = jpiglo-i1+1
560         ENDIF
561      ELSE IF (isens ==2) THEN
562         IF (pty == 2) THEN ! T, U points
563            agrif_external_switch_index = jpjglo-2-(i1 -jpjglo)
564         ELSE ! V, F points
565            agrif_external_switch_index = jpjglo-3-(i1 -jpjglo)
566         ENDIF
567      ENDIF
568
569   END FUNCTION agrif_external_switch_index
570
571   SUBROUTINE correct_field(tab2d,i1,i2,j1,j2)
572      USE dom_oce
573      INTEGER :: i1,i2,j1,j2
574      REAL,DIMENSION(i1:i2,j1:j2) :: tab2d
575
576      INTEGER :: i,j
577      REAL,DIMENSION(i1:i2,j1:j2) :: tab2dtemp
578
579      tab2dtemp = tab2d
580
581      DO j=j1,j2
582         DO i=i1,i2
583        tab2d(i,j)=tab2dtemp(i2-(i-i1),j2-(j-j1))
584         END DO
585      ENDDO
586
587   END SUBROUTINE correct_field
588
589   SUBROUTINE agrif_init_lonlat()
590      USE agrif_profiles
591      USE agrif_util
592      USE dom_oce
593   
594      LOGICAL  :: l_deg 
595      EXTERNAL :: init_glamt, init_glamu, init_glamv, init_glamf
596      EXTERNAL :: init_gphit, init_gphiu, init_gphiv, init_gphif
597      REAL,EXTERNAL :: longitude_linear_interp
598
599      l_deg = .TRUE.
600      IF  ( Agrif_Parent(jphgr_msh)==2 & 
601      & .OR.Agrif_Parent(jphgr_msh)==3 & 
602      & .OR.Agrif_Parent(jphgr_msh)==5 ) l_deg = .FALSE.
603
604      IF ( l_deg ) THEN
605         CALL Agrif_Set_external_linear_interp(longitude_linear_interp)
606      ENDIF
607
608      CALL Agrif_Init_variable(glamt_id, procname = init_glamt)
609      CALL Agrif_Init_variable(glamu_id, procname = init_glamu)
610      CALL Agrif_Init_variable(glamv_id, procname = init_glamv)
611      CALL Agrif_Init_variable(glamf_id, procname = init_glamf)
612      CALL Agrif_UnSet_external_linear_interp()
613
614      CALL Agrif_Init_variable(gphit_id, procname = init_gphit)
615      CALL Agrif_Init_variable(gphiu_id, procname = init_gphiu)
616      CALL Agrif_Init_variable(gphiv_id, procname = init_gphiv)
617      CALL Agrif_Init_variable(gphif_id, procname = init_gphif)
618
619   END SUBROUTINE agrif_init_lonlat
620
621   REAL FUNCTION longitude_linear_interp(x1,x2,coeff)
622      REAL :: x1, x2, coeff
623      REAL :: val_interp
624
625      IF( (x1*x2 <= -50*50) ) THEN
626         IF( x1 < 0 ) THEN
627            val_interp = coeff *(x1+360.) + (1.-coeff) *x2
628         ELSE
629            val_interp = coeff *x1 + (1.-coeff) *(x2+360.)
630         ENDIF
631         IF ((val_interp) >=180.) val_interp = val_interp - 360.
632      ELSE
633      val_interp = coeff * x1 + (1.-coeff) * x2
634      ENDIF
635      longitude_linear_interp = val_interp
636
637   END FUNCTION longitude_linear_interp
638
639   SUBROUTINE agrif_init_scales()
640      USE agrif_profiles
641      USE agrif_util
642      USE dom_oce
643      USE lbclnk
644      LOGICAL :: ln_perio
645      INTEGER jpi,jpj
646
647      EXTERNAL :: init_e1t, init_e1u, init_e1v, init_e1f
648      EXTERNAL :: init_e2t, init_e2u, init_e2v, init_e2f
649
650      ln_perio=.FALSE.
651      if( jperio ==1 .OR. jperio==2 .OR. jperio==4) ln_perio=.TRUE.
652
653      CALL Agrif_Init_variable(e1t_id, procname = init_e1t)
654      CALL Agrif_Init_variable(e1u_id, procname = init_e1u)
655      CALL Agrif_Init_variable(e1v_id, procname = init_e1v)
656      CALL Agrif_Init_variable(e1f_id, procname = init_e1f)
657
658      CALL Agrif_Init_variable(e2t_id, procname = init_e2t)
659      CALL Agrif_Init_variable(e2u_id, procname = init_e2u)
660      CALL Agrif_Init_variable(e2v_id, procname = init_e2v)
661      CALL Agrif_Init_variable(e2f_id, procname = init_e2f)
662
663      CALL lbc_lnk( 'e1t', e1t, 'T', 1._wp)
664      CALL lbc_lnk( 'e2t', e2t, 'T', 1._wp)
665      CALL lbc_lnk( 'e1u', e1u, 'U', 1._wp)
666      CALL lbc_lnk( 'e2u', e2u, 'U', 1._wp)
667      CALL lbc_lnk( 'e1v', e1v, 'V', 1._wp)
668      CALL lbc_lnk( 'e2v', e2v, 'V', 1._wp)
669      CALL lbc_lnk( 'e1f', e1f, 'F', 1._wp)
670      CALL lbc_lnk( 'e2f', e2f, 'F', 1._wp)
671
672   END SUBROUTINE agrif_init_scales
673
674   SUBROUTINE init_glamt( ptab, i1, i2, j1, j2,  before, nb,ndir)
675      USE dom_oce
676      !!----------------------------------------------------------------------
677      !!                  ***  ROUTINE interpsshn  ***
678      !!----------------------------------------------------------------------
679      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
680      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
681      LOGICAL                         , INTENT(in   ) ::   before
682      INTEGER                         , INTENT(in   ) ::   nb , ndir
683
684      !
685      !!----------------------------------------------------------------------
686      !
687      IF( before) THEN
688         ptab(i1:i2,j1:j2) = glamt(i1:i2,j1:j2)
689      ELSE
690          glamt(i1:i2,j1:j2) = ptab(i1:i2,j1:j2)
691      ENDIF
692      !
693   END SUBROUTINE init_glamt
694
695   SUBROUTINE init_glamu( ptab, i1, i2, j1, j2, before, nb,ndir)
696      USE dom_oce
697      !!----------------------------------------------------------------------
698      !!                  ***  ROUTINE interpsshn  ***
699      !!----------------------------------------------------------------------
700      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
701      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
702      LOGICAL                         , INTENT(in   ) ::   before
703      INTEGER                         , INTENT(in   ) ::   nb , ndir
704      LOGICAL  ::   western_side, eastern_side,northern_side,southern_side
705      !
706      !!----------------------------------------------------------------------
707      !
708      IF( before) THEN
709         ptab(i1:i2,j1:j2) = glamu(i1:i2,j1:j2)
710      ELSE
711          glamu(i1:i2,j1:j2) = ptab(i1:i2,j1:j2)
712      ENDIF
713      !
714   END SUBROUTINE init_glamu
715
716   SUBROUTINE init_glamv( ptab, i1, i2, j1, j2, before, nb,ndir)
717      USE dom_oce
718      !!----------------------------------------------------------------------
719      !!                  ***  ROUTINE interpsshn  ***
720      !!----------------------------------------------------------------------
721      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
722      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
723      LOGICAL                         , INTENT(in   ) ::   before
724      INTEGER                         , INTENT(in   ) ::   nb , ndir
725      !
726      !!----------------------------------------------------------------------
727      !
728      IF( before) THEN
729         ptab(i1:i2,j1:j2) = glamv(i1:i2,j1:j2)
730      ELSE
731          glamv(i1:i2,j1:j2) = ptab(i1:i2,j1:j2)
732      ENDIF
733      !
734   END SUBROUTINE init_glamv
735
736   SUBROUTINE init_glamf( ptab, i1, i2, j1, j2,  before, nb,ndir)
737      USE dom_oce
738      !!----------------------------------------------------------------------
739      !!                  ***  ROUTINE init_glamf  ***
740      !!----------------------------------------------------------------------
741      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
742      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
743      LOGICAL                         , INTENT(in   ) ::   before
744      INTEGER                         , INTENT(in   ) ::   nb , ndir
745      !
746      !!----------------------------------------------------------------------
747      !
748      IF( before) THEN
749         ptab(i1:i2,j1:j2) = glamf(i1:i2,j1:j2)
750      ELSE
751          glamf(i1:i2,j1:j2) = ptab(i1:i2,j1:j2)
752      ENDIF
753      !
754   END SUBROUTINE init_glamf
755
756   SUBROUTINE init_gphit( ptab, i1, i2, j1, j2, before, nb,ndir)
757      USE dom_oce
758      !!----------------------------------------------------------------------
759      !!                  ***  ROUTINE init_gphit  ***
760      !!----------------------------------------------------------------------
761      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
762      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
763      LOGICAL                         , INTENT(in   ) ::   before
764      INTEGER                         , INTENT(in   ) ::   nb , ndir
765      !
766      !!----------------------------------------------------------------------
767      !
768      IF( before) THEN
769         ptab(i1:i2,j1:j2) = gphit(i1:i2,j1:j2)
770      ELSE
771         gphit(i1:i2,j1:j2)=ptab
772      ENDIF
773      !
774   END SUBROUTINE init_gphit
775
776   SUBROUTINE init_gphiu( ptab, i1, i2, j1, j2, before, nb,ndir)
777      USE dom_oce
778      !!----------------------------------------------------------------------
779      !!                  ***  ROUTINE init_gphiu  ***
780      !!----------------------------------------------------------------------
781      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
782      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
783      LOGICAL                         , INTENT(in   ) ::   before
784      INTEGER                         , INTENT(in   ) ::   nb , ndir
785      !
786      !!----------------------------------------------------------------------
787      !
788      IF( before) THEN
789         ptab(i1:i2,j1:j2) = gphiu(i1:i2,j1:j2)
790      ELSE
791         gphiu(i1:i2,j1:j2)=ptab
792      ENDIF
793      !
794   END SUBROUTINE init_gphiu
795
796   SUBROUTINE init_gphiv( ptab, i1, i2, j1, j2, before, nb,ndir)
797      USE dom_oce
798      !!----------------------------------------------------------------------
799      !!                  ***  ROUTINE init_gphiv  ***
800      !!----------------------------------------------------------------------
801      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
802      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
803      LOGICAL                         , INTENT(in   ) ::   before
804      INTEGER                         , INTENT(in   ) ::   nb , ndir
805      !
806      !!----------------------------------------------------------------------
807
808      IF( before) THEN
809         ptab(i1:i2,j1:j2) = gphiv(i1:i2,j1:j2)
810      ELSE
811         gphiv(i1:i2,j1:j2)=ptab
812      ENDIF
813      !
814   END SUBROUTINE init_gphiv
815
816
817   SUBROUTINE init_gphif( ptab, i1, i2, j1, j2, before, nb,ndir)
818      USE dom_oce
819      !!----------------------------------------------------------------------
820      !!                  ***  ROUTINE init_gphif  ***
821      !!----------------------------------------------------------------------
822      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
823      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
824      LOGICAL                         , INTENT(in   ) ::   before
825      INTEGER                         , INTENT(in   ) ::   nb , ndir
826      !
827      !!----------------------------------------------------------------------
828      !
829      IF( before) THEN
830         ptab(i1:i2,j1:j2) = gphif(i1:i2,j1:j2)
831      ELSE
832         gphif(i1:i2,j1:j2)=ptab
833      ENDIF
834      !
835   END SUBROUTINE init_gphif
836
837
838   SUBROUTINE init_e1t( ptab, i1, i2, j1, j2, before, nb,ndir)
839      USE dom_oce
840      USE agrif_parameters
841      !!----------------------------------------------------------------------
842      !!                  ***  ROUTINE init_e1t  ***
843      !!----------------------------------------------------------------------
844      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
845      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
846      LOGICAL                         , INTENT(in   ) ::   before
847      INTEGER                         , INTENT(in   ) ::   nb , ndir
848      !
849      !!----------------------------------------------------------------------
850      !
851      INTEGER :: jj
852
853      IF( before) THEN
854        ! May need to extend at south boundary
855          IF (j1<1) THEN
856            IF (.NOT.agrif_child(lk_south)) THEN
857              IF ((nbondj == -1).OR.(nbondj == 2)) THEN
858                DO jj=1,j2
859                  ptab(i1:i2,jj)=e1t(i1:i2,jj)
860                ENDDO
861                DO jj=j1,0
862                  ptab(i1:i2,jj)=e1t(i1:i2,1)
863                ENDDO
864              ENDIF
865            ELSE
866              stop "OUT OF BOUNDS"
867            ENDIF
868          ELSE
869             ptab(i1:i2,j1:j2) = e1t(i1:i2,j1:j2)
870          ENDIF
871      ELSE
872         e1t(i1:i2,j1:j2)=ptab/Agrif_rhoy()
873      ENDIF
874      !
875   END SUBROUTINE init_e1t
876
877   SUBROUTINE init_e1u( ptab, i1, i2, j1, j2, before, nb,ndir)
878      USE dom_oce
879      USE agrif_parameters
880      !!----------------------------------------------------------------------
881      !!                  ***  ROUTINE init_e1u  ***
882      !!----------------------------------------------------------------------
883      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
884      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
885      LOGICAL                         , INTENT(in   ) ::   before
886      INTEGER                         , INTENT(in   ) ::   nb , ndir
887      !
888      !!----------------------------------------------------------------------
889      !
890      INTEGER :: jj
891
892      IF( before) THEN
893          IF (j1<1) THEN
894            IF (.NOT.agrif_child(lk_south)) THEN
895              IF ((nbondj == -1).OR.(nbondj == 2)) THEN
896                DO jj=1,j2
897                  ptab(i1:i2,jj)=e1u(i1:i2,jj)
898                ENDDO
899                DO jj=j1,0
900                  ptab(i1:i2,jj)=e1u(i1:i2,1)
901                ENDDO
902              ENDIF
903            ELSE
904              stop "OUT OF BOUNDS"
905            ENDIF
906          ELSE
907             ptab(i1:i2,j1:j2) = e1u(i1:i2,j1:j2)
908          ENDIF
909      ELSE
910         e1u(i1:i2,j1:j2)=ptab/Agrif_rhoy()
911      ENDIF
912      !
913   END SUBROUTINE init_e1u
914
915   SUBROUTINE init_e1v( ptab, i1, i2, j1, j2, before, nb,ndir)
916      USE dom_oce
917      !!----------------------------------------------------------------------
918      !!                  ***  ROUTINE init_e1v  ***
919      !!----------------------------------------------------------------------
920      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
921      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
922      LOGICAL                         , INTENT(in   ) ::   before
923      INTEGER                         , INTENT(in   ) ::   nb , ndir
924      !
925      !!----------------------------------------------------------------------
926      !
927      IF( before) THEN
928         ptab(i1:i2,j1:j2) = e1v(i1:i2,j1:j2)
929      ELSE
930         e1v(i1:i2,j1:j2)=ptab/Agrif_rhoy()
931      ENDIF
932      !
933   END SUBROUTINE init_e1v
934
935   SUBROUTINE init_e1f( ptab, i1, i2, j1, j2, before, nb,ndir)
936      USE dom_oce
937      !!----------------------------------------------------------------------
938      !!                  ***  ROUTINE init_e1f  ***
939      !!----------------------------------------------------------------------
940      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
941      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
942      LOGICAL                         , INTENT(in   ) ::   before
943      INTEGER                         , INTENT(in   ) ::   nb , ndir
944      !
945      !!----------------------------------------------------------------------
946      !
947      IF( before) THEN
948         ptab(i1:i2,j1:j2) = e1f(i1:i2,j1:j2)
949      ELSE
950         e1f(i1:i2,j1:j2)=ptab/Agrif_rhoy()
951      ENDIF
952      !
953   END SUBROUTINE init_e1f
954
955   SUBROUTINE init_e2t( ptab, i1, i2, j1, j2, before, nb,ndir)
956      USE dom_oce
957      USE agrif_parameters
958      !!----------------------------------------------------------------------
959      !!                  ***  ROUTINE init_e2t  ***
960      !!----------------------------------------------------------------------
961      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
962      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
963      LOGICAL                         , INTENT(in   ) ::   before
964      INTEGER                         , INTENT(in   ) ::   nb , ndir
965      !
966      !!----------------------------------------------------------------------
967      !
968      INTEGER :: jj
969
970      IF( before) THEN
971          IF (j1<1) THEN
972            IF (.NOT.agrif_child(lk_south)) THEN
973              IF ((nbondj == -1).OR.(nbondj == 2)) THEN
974                DO jj=1,j2
975                  ptab(i1:i2,jj)=e2t(i1:i2,jj)
976                ENDDO
977                DO jj=j1,0
978                  ptab(i1:i2,jj)=e2t(i1:i2,1)
979                ENDDO
980              ENDIF
981            ELSE
982              stop "OUT OF BOUNDS"
983            ENDIF
984          ELSE
985             ptab(i1:i2,j1:j2) = e2t(i1:i2,j1:j2)
986          ENDIF
987      ELSE
988         e2t(i1:i2,j1:j2)=ptab/Agrif_rhoy()
989      ENDIF
990      !
991   END SUBROUTINE init_e2t
992
993   SUBROUTINE init_e2u( ptab, i1, i2, j1, j2, before, nb,ndir)
994   USE dom_oce
995   USE agrif_parameters
996      !!----------------------------------------------------------------------
997      !!                  ***  ROUTINE interpsshn  ***
998      !!----------------------------------------------------------------------
999      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
1000      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
1001      LOGICAL                         , INTENT(in   ) ::   before
1002      INTEGER                         , INTENT(in   ) ::   nb , ndir
1003      !
1004      !!----------------------------------------------------------------------
1005      !
1006      INTEGER :: jj
1007
1008      IF( before) THEN
1009          IF (j1<1) THEN
1010            IF (.NOT.agrif_child(lk_south)) THEN
1011              IF ((nbondj == -1).OR.(nbondj == 2)) THEN
1012                DO jj=1,j2
1013                  ptab(i1:i2,jj)=e2u(i1:i2,jj)
1014                ENDDO
1015                DO jj=j1,0
1016                  ptab(i1:i2,jj)=e2u(i1:i2,1)
1017                ENDDO
1018              ENDIF
1019            ELSE
1020              stop "OUT OF BOUNDS"
1021            ENDIF
1022          ELSE
1023             ptab(i1:i2,j1:j2) = e2u(i1:i2,j1:j2)
1024          ENDIF
1025      ELSE
1026         e2u(i1:i2,j1:j2)=ptab/Agrif_rhoy()
1027      ENDIF
1028      !
1029   END SUBROUTINE init_e2u
1030
1031   SUBROUTINE init_e2v( ptab, i1, i2, j1, j2, before, nb,ndir)
1032      USE dom_oce
1033      !!----------------------------------------------------------------------
1034      !!                  ***  ROUTINE interpsshn  ***
1035      !!----------------------------------------------------------------------
1036      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
1037      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
1038      LOGICAL                         , INTENT(in   ) ::   before
1039      INTEGER                         , INTENT(in   ) ::   nb , ndir
1040      !
1041      !!----------------------------------------------------------------------
1042      IF( before) THEN
1043         ptab(i1:i2,j1:j2) = e2v(i1:i2,j1:j2)
1044      ELSE
1045         e2v(i1:i2,j1:j2)=ptab/Agrif_rhoy()
1046      ENDIF
1047      !
1048   END SUBROUTINE init_e2v
1049
1050   SUBROUTINE init_e2f( ptab, i1, i2, j1, j2, before, nb,ndir)
1051   USE dom_oce
1052      !!----------------------------------------------------------------------
1053      !!                  ***  ROUTINE interpsshn  ***
1054      !!----------------------------------------------------------------------
1055      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2
1056      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab
1057      LOGICAL                         , INTENT(in   ) ::   before
1058      INTEGER                         , INTENT(in   ) ::   nb , ndir
1059      !
1060      !!----------------------------------------------------------------------
1061      !
1062      IF( before) THEN
1063         ptab(i1:i2,j1:j2) = e2f(i1:i2,j1:j2)
1064      ELSE
1065         e2f(i1:i2,j1:j2)=ptab/Agrif_rhoy()
1066      ENDIF
1067      !
1068   END SUBROUTINE init_e2f
1069
1070
1071   SUBROUTINE agrif_nemo_init
1072      USE agrif_parameters
1073      USE dom_oce
1074      USE in_out_manager
1075      USE lib_mpp
1076      !!
1077      IMPLICIT NONE
1078
1079      INTEGER ::   ios
1080
1081      NAMELIST/namagrif/ nn_cln_update,ln_spc_dyn,rn_sponge_tra,rn_sponge_dyn,ln_chk_bathy,npt_connect,   &
1082      &  npt_copy
1083
1084  !    REWIND( numnam_ref )              ! Namelist namagrif in reference namelist : nesting parameters
1085      READ  ( numnam_ref, namagrif, IOSTAT = ios, ERR = 901 )
1086901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namagrif in reference namelist')
1087
1088  !    REWIND( numnam_cfg )              ! Namelist namzgr in configuration namelist : nesting parameters
1089      READ  ( numnam_cfg, namagrif, IOSTAT = ios, ERR = 902 )
1090902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namagrif in configuration namelist')
1091      IF(lwm) WRITE ( numond, namagrif )
1092
1093      IF(lwp) THEN                     ! Control print
1094         WRITE(numout,*)
1095         WRITE(numout,*) 'agrif_nemo_init : nesting'
1096         WRITE(numout,*) '~~~~~~~'
1097         WRITE(numout,*) '   Namelist namagrif : set nesting parameters'
1098         WRITE(numout,*) '      npt_copy             = ', npt_copy
1099         WRITE(numout,*) '      npt_connect          = ', npt_connect
1100      ENDIF
1101
1102   ! Set the number of ghost cells according to periodicity
1103
1104      nbghostcells_x = nbghostcells
1105      nbghostcells_y_s = nbghostcells
1106      nbghostcells_y_n = nbghostcells
1107
1108      IF ((jperio == 1).OR.(jperio == 4)) THEN
1109        nbghostcells_x = 0
1110      ENDIF
1111      IF (jperio == 4) THEN
1112        nbghostcells_y_s = 0
1113      ENDIF
1114
1115      IF (.not.agrif_root()) THEN
1116         lk_west  = .NOT. ( Agrif_Ix() == 1 )
1117         lk_east  = .NOT. ( Agrif_Ix() + nbcellsx/AGRIF_Irhox() == Agrif_Parent(Ni0glo) + 1 )
1118         lk_south = .NOT. ( Agrif_Iy() == 1 )
1119         lk_north = .NOT. ( Agrif_Iy() + nbcellsy/AGRIF_Irhoy() == Agrif_Parent(Nj0glo) + 1)
1120         IF (.NOT.lk_south) THEN
1121            nbghostcells_y_s = 0
1122         ENDIF
1123         IF (.NOT.lk_north) THEN
1124            nbghostcells_y_n = 0
1125         ENDIF
1126         IF(lwp) THEN                     ! Control print
1127            WRITE(numout,*)
1128            WRITE(numout,*) 'nbghostcells_y_s', nbghostcells_y_s
1129            WRITE(numout,*) 'nbghostcells_y_n', nbghostcells_y_n
1130            WRITE(numout,*) 'nbghostcells_x', nbghostcells_x
1131            WRITE(numout,*) 'lk_west', lk_west
1132            WRITE(numout,*) 'lk_east', lk_east
1133            WRITE(numout,*) 'lk_south', lk_south
1134            WRITE(numout,*) 'lk_north', lk_north
1135         ENDIF
1136      ELSE ! root grid
1137         nbghostcells_x   = 0 
1138         nbghostcells_y_s = 0 
1139         nbghostcells_y_n = 0 
1140      ENDIF
1141
1142   END SUBROUTINE agrif_nemo_init
1143
1144   SUBROUTINE Agrif_detect( kg, ksizex )
1145      !!----------------------------------------------------------------------
1146      !!                      *** ROUTINE Agrif_detect ***
1147      !!----------------------------------------------------------------------
1148      INTEGER, DIMENSION(2) :: ksizex
1149      INTEGER, DIMENSION(ksizex(1),ksizex(2)) :: kg
1150      !!----------------------------------------------------------------------
1151      !
1152      RETURN
1153      !
1154   END SUBROUTINE Agrif_detect
1155
1156   SUBROUTINE agrif_before_regridding
1157   END SUBROUTINE agrif_before_regridding
1158
1159# if defined  key_mpp_mpi
1160   SUBROUTINE Agrif_InvLoc( indloc, nprocloc, i, indglob )
1161         !!----------------------------------------------------------------------
1162         !!                     *** ROUTINE Agrif_InvLoc ***
1163         !!----------------------------------------------------------------------
1164      USE dom_oce
1165      !!
1166      IMPLICIT NONE
1167      !
1168      INTEGER :: indglob, indloc, nprocloc, i
1169         !!----------------------------------------------------------------------
1170      !
1171      SELECT CASE( i )
1172      CASE(1)   ;   indglob = mig(indloc)
1173      CASE(2)   ;   indglob = mjg(indloc)
1174      CASE DEFAULT
1175         indglob = indloc
1176      END SELECT
1177      !
1178   END SUBROUTINE Agrif_InvLoc
1179
1180   SUBROUTINE Agrif_get_proc_info( imin, imax, jmin, jmax )
1181      !!----------------------------------------------------------------------
1182      !!                 *** ROUTINE Agrif_get_proc_info ***
1183      !!----------------------------------------------------------------------
1184      USE par_oce
1185      USE dom_oce
1186      !!
1187      IMPLICIT NONE
1188      !
1189      INTEGER, INTENT(out) :: imin, imax
1190      INTEGER, INTENT(out) :: jmin, jmax
1191         !!----------------------------------------------------------------------
1192      !
1193      imin = nimppt(Agrif_Procrank+1)  ! ?????
1194      jmin = njmppt(Agrif_Procrank+1)  ! ?????
1195      imax = imin + jpi - 1
1196      jmax = jmin + jpj - 1
1197      !
1198   END SUBROUTINE Agrif_get_proc_info
1199
1200   SUBROUTINE Agrif_estimate_parallel_cost(imin, imax,jmin, jmax, nbprocs, grid_cost)
1201      !!----------------------------------------------------------------------
1202      !!                 *** ROUTINE Agrif_estimate_parallel_cost ***
1203      !!----------------------------------------------------------------------
1204      USE par_oce
1205      !!
1206      IMPLICIT NONE
1207      !
1208      INTEGER,  INTENT(in)  :: imin, imax
1209      INTEGER,  INTENT(in)  :: jmin, jmax
1210      INTEGER,  INTENT(in)  :: nbprocs
1211      REAL(wp), INTENT(out) :: grid_cost
1212         !!----------------------------------------------------------------------
1213      !
1214      grid_cost = REAL(imax-imin+1,wp)*REAL(jmax-jmin+1,wp) / REAL(nbprocs,wp)
1215      !
1216   END SUBROUTINE Agrif_estimate_parallel_cost
1217# endif
1218#endif
Note: See TracBrowser for help on using the repository browser.