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

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

#2638, use i refinement factor in place of j to define e1x over child grids

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