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 @ 14623

Last change on this file since 14623 was 14623, checked in by ldebreu, 3 years ago

AGFdomcfg: 1) Update DOMAINcfg to be compliant with the removal of halo cells 2) Update most of the LBC ... subroutines to a recent NEMO 4 version #2638

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