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