1 |
#define CPL_DEBUG |
2 |
#define FIX_FOR_EDGE_WINDS |
3 |
#include "PACKAGES_CONFIG.h" |
4 |
#include "CPP_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: CPL_MPMICE |
8 |
C !INTERFACE: |
9 |
SUBROUTINE CPL_MPMICE( myTime, myIter, myThid ) |
10 |
|
11 |
C !DESCRIPTION: \bv |
12 |
C *================================================================== |
13 |
C | SUBROUTINE cpl_mpmice |
14 |
C | o Couple MITgcm ocean model with MPMice sea ice model |
15 |
C *================================================================== |
16 |
C \ev |
17 |
|
18 |
C !USES: |
19 |
IMPLICIT NONE |
20 |
C == Global variables == |
21 |
#include "SIZE.h" |
22 |
#include "EEPARAMS.h" |
23 |
#include "PARAMS.h" |
24 |
#include "DYNVARS.h" |
25 |
#include "GRID.h" |
26 |
#ifdef ALLOW_EXF |
27 |
# include "EXF_OPTIONS.h" |
28 |
# include "EXF_FIELDS.h" |
29 |
#endif |
30 |
#ifdef ALLOW_SEAICE |
31 |
# include "SEAICE_OPTIONS.h" |
32 |
# include "SEAICE.h" |
33 |
#endif |
34 |
|
35 |
LOGICAL DIFFERENT_MULTIPLE |
36 |
EXTERNAL DIFFERENT_MULTIPLE |
37 |
|
38 |
C !LOCAL VARIABLES: |
39 |
C mytime - time counter for this thread (seconds) |
40 |
C myiter - iteration counter for this thread |
41 |
C mythid - thread number for this instance of the routine. |
42 |
_RL mytime |
43 |
INTEGER myiter, mythid |
44 |
CEOP |
45 |
|
46 |
#ifdef ALLOW_CPL_MPMICE |
47 |
# include "EESUPPORT.h" |
48 |
# include "CPL_MPMICE.h" |
49 |
COMMON /CPL_MPI_ID/ |
50 |
& myworldid, local_ocean_leader, local_ice_leader |
51 |
integer myworldid, local_ocean_leader, local_ice_leader |
52 |
integer mpistatus(MPI_STATUS_SIZE), mpierr |
53 |
integer xfer_gridsize(2) |
54 |
integer i, j, bi, bj, buffsize, idx |
55 |
Real*8 xfer_scalar |
56 |
Real*8 xfer_array(Nx,Ny) |
57 |
Real*8 xfer_bc_tracer(2*(Nx+Ny)-4) |
58 |
Real*8 xfer_bc_veloc(2*(Nx+Ny)-6) |
59 |
_RL local(1:sNx,1:sNy,nSx,nSy) |
60 |
_RL ScatArray(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
61 |
|
62 |
IF( myTime .EQ. startTime ) THEN |
63 |
|
64 |
C Send deltatimestep |
65 |
_BEGIN_MASTER( myThid ) |
66 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
67 |
xfer_scalar = deltat |
68 |
buffsize = 1 |
69 |
CALL MPI_SEND(xfer_scalar,buffsize,MPI_DOUBLE_PRECISION, |
70 |
& local_ice_leader,TimeIntervalTag,MPI_COMM_WORLD,mpierr) |
71 |
#ifdef CPL_DEBUG |
72 |
print*,'MITgcm send TimeInterval', xfer_scalar |
73 |
#endif |
74 |
ENDIF |
75 |
_END_MASTER( myThid ) |
76 |
|
77 |
C Send grid dimensions (Nx,Ny) |
78 |
_BEGIN_MASTER( myThid ) |
79 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
80 |
xfer_gridsize(1)=Nx |
81 |
xfer_gridsize(2)=Ny |
82 |
buffsize = 2 |
83 |
CALL MPI_SEND(xfer_gridsize,buffsize,MPI_INTEGER, |
84 |
& local_ice_leader,OceanGridsizeTag,MPI_COMM_WORLD,mpierr) |
85 |
#ifdef CPL_DEBUG |
86 |
print*,'MITgcm send OceanGridsize', xfer_gridsize |
87 |
#endif |
88 |
ENDIF |
89 |
_END_MASTER( myThid ) |
90 |
|
91 |
C Send ice area |
92 |
DO bj=1,nSy |
93 |
DO bi=1,nSx |
94 |
DO j=1,sNy |
95 |
DO i=1,sNx |
96 |
local(i,j,bi,bj) = AREA(i,j,1,bi,bj) |
97 |
ENDDO |
98 |
ENDDO |
99 |
ENDDO |
100 |
ENDDO |
101 |
CALL GATHER_2D( xfer_array, local, myThid ) |
102 |
_BEGIN_MASTER( myThid ) |
103 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
104 |
buffsize = Nx*Ny |
105 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
106 |
& local_ice_leader,AreaTag,MPI_COMM_WORLD,mpierr) |
107 |
ENDIF |
108 |
_END_MASTER( myThid ) |
109 |
#ifdef CPL_DEBUG |
110 |
CALL PLOT_FIELD_XYRL( AREA, 'AREA', myIter, myThid ) |
111 |
#endif |
112 |
|
113 |
C Send ice thickness |
114 |
DO bj=1,nSy |
115 |
DO bi=1,nSx |
116 |
DO j=1,sNy |
117 |
DO i=1,sNx |
118 |
local(i,j,bi,bj) = HEFF(i,j,1,bi,bj) |
119 |
ENDDO |
120 |
ENDDO |
121 |
ENDDO |
122 |
ENDDO |
123 |
CALL GATHER_2D( xfer_array, local, myThid ) |
124 |
_BEGIN_MASTER( myThid ) |
125 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
126 |
buffsize = Nx*Ny |
127 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
128 |
& local_ice_leader,HeffTag,MPI_COMM_WORLD,mpierr) |
129 |
ENDIF |
130 |
_END_MASTER( myThid ) |
131 |
#ifdef CPL_DEBUG |
132 |
CALL PLOT_FIELD_XYRL( HEFF, 'HEFF', myIter, myThid ) |
133 |
#endif |
134 |
|
135 |
C Send ice salinity |
136 |
DO bj=1,nSy |
137 |
DO bi=1,nSx |
138 |
DO j=1,sNy |
139 |
DO i=1,sNx |
140 |
local(i,j,bi,bj) = HSALT(i,j,bi,bj) |
141 |
ENDDO |
142 |
ENDDO |
143 |
ENDDO |
144 |
ENDDO |
145 |
CALL GATHER_2D( xfer_array, local, myThid ) |
146 |
_BEGIN_MASTER( myThid ) |
147 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
148 |
buffsize = Nx*Ny |
149 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
150 |
& local_ice_leader,HsaltTag,MPI_COMM_WORLD,mpierr) |
151 |
ENDIF |
152 |
_END_MASTER( myThid ) |
153 |
#ifdef CPL_DEBUG |
154 |
CALL PLOT_FIELD_XYRL( HSALT, 'HSALT', myIter, myThid ) |
155 |
#endif |
156 |
|
157 |
C Send snow thickness |
158 |
DO bj=1,nSy |
159 |
DO bi=1,nSx |
160 |
DO j=1,sNy |
161 |
DO i=1,sNx |
162 |
local(i,j,bi,bj) = HSNOW(i,j,bi,bj) |
163 |
ENDDO |
164 |
ENDDO |
165 |
ENDDO |
166 |
ENDDO |
167 |
CALL GATHER_2D( xfer_array, local, myThid ) |
168 |
_BEGIN_MASTER( myThid ) |
169 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
170 |
buffsize = Nx*Ny |
171 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
172 |
& local_ice_leader,HsnowTag,MPI_COMM_WORLD,mpierr) |
173 |
ENDIF |
174 |
_END_MASTER( myThid ) |
175 |
#ifdef CPL_DEBUG |
176 |
CALL PLOT_FIELD_XYRL( HSNOW, 'HSNOW', myIter, myThid ) |
177 |
#endif |
178 |
|
179 |
ENDIF ! ( myTime .EQ. startTime ) |
180 |
|
181 |
C Send ocean model time |
182 |
_BEGIN_MASTER( myThid ) |
183 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
184 |
xfer_scalar = myTime |
185 |
buffsize = 1 |
186 |
CALL MPI_SEND(xfer_scalar,buffsize,MPI_DOUBLE_PRECISION, |
187 |
& local_ice_leader,OceanTimeTag,MPI_COMM_WORLD,mpierr) |
188 |
#ifdef CPL_DEBUG |
189 |
print*,'MITgcm send OceanTime', xfer_scalar |
190 |
#endif |
191 |
ENDIF |
192 |
_END_MASTER( myThid ) |
193 |
|
194 |
C Send boundary ice area |
195 |
DO bj=1,nSy |
196 |
DO bi=1,nSx |
197 |
DO j=1,sNy |
198 |
DO i=1,sNx |
199 |
local(i,j,bi,bj) = AREA(i,j,1,bi,bj) |
200 |
ENDDO |
201 |
ENDDO |
202 |
ENDDO |
203 |
ENDDO |
204 |
CALL GATHER_2D( xfer_array, local, myThid ) |
205 |
idx = 0 |
206 |
DO i = 1, Nx |
207 |
idx = idx + 1 |
208 |
xfer_bc_tracer(idx) = xfer_array(i,1) |
209 |
ENDDO |
210 |
DO j = 2, Ny |
211 |
idx = idx + 1 |
212 |
xfer_bc_tracer(idx) = xfer_array(Nx,j) |
213 |
ENDDO |
214 |
DO i = (Nx-1), -1, 1 |
215 |
idx = idx + 1 |
216 |
xfer_bc_tracer(idx) = xfer_array(i,Ny) |
217 |
ENDDO |
218 |
DO j = (Ny-1), -1, 2 |
219 |
idx = idx + 1 |
220 |
xfer_bc_tracer(idx) = xfer_array(1,j) |
221 |
ENDDO |
222 |
_BEGIN_MASTER( myThid ) |
223 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
224 |
buffsize = 2*(Nx+Ny)-4 |
225 |
print*,'MITgcm is about to send AreaBcTag',buffsize |
226 |
CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION, |
227 |
& local_ice_leader,AreaBcTag,MPI_COMM_WORLD,mpierr) |
228 |
print*,'MITgcm has sent AreaBcTag',buffsize |
229 |
ENDIF |
230 |
_END_MASTER( myThid ) |
231 |
|
232 |
C Send boundary ice thickness |
233 |
DO bj=1,nSy |
234 |
DO bi=1,nSx |
235 |
DO j=1,sNy |
236 |
DO i=1,sNx |
237 |
local(i,j,bi,bj) = HEFF(i,j,1,bi,bj) |
238 |
ENDDO |
239 |
ENDDO |
240 |
ENDDO |
241 |
ENDDO |
242 |
CALL GATHER_2D( xfer_array, local, myThid ) |
243 |
idx = 0 |
244 |
DO i = 1, Nx |
245 |
idx = idx + 1 |
246 |
xfer_bc_tracer(idx) = xfer_array(i,1) |
247 |
ENDDO |
248 |
DO j = 2, Ny |
249 |
idx = idx + 1 |
250 |
xfer_bc_tracer(idx) = xfer_array(Nx,j) |
251 |
ENDDO |
252 |
DO i = (Nx-1), -1, 1 |
253 |
idx = idx + 1 |
254 |
xfer_bc_tracer(idx) = xfer_array(i,Ny) |
255 |
ENDDO |
256 |
DO j = (Ny-1), -1, 2 |
257 |
idx = idx + 1 |
258 |
xfer_bc_tracer(idx) = xfer_array(1,j) |
259 |
ENDDO |
260 |
_BEGIN_MASTER( myThid ) |
261 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
262 |
buffsize = 2*(Nx+Ny)-4 |
263 |
CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION, |
264 |
& local_ice_leader,HeffBcTag,MPI_COMM_WORLD,mpierr) |
265 |
ENDIF |
266 |
_END_MASTER( myThid ) |
267 |
|
268 |
C Send boundary ice salinity |
269 |
DO bj=1,nSy |
270 |
DO bi=1,nSx |
271 |
DO j=1,sNy |
272 |
DO i=1,sNx |
273 |
local(i,j,bi,bj) = HSALT(i,j,bi,bj) |
274 |
ENDDO |
275 |
ENDDO |
276 |
ENDDO |
277 |
ENDDO |
278 |
CALL GATHER_2D( xfer_array, local, myThid ) |
279 |
idx = 0 |
280 |
DO i = 1, Nx |
281 |
idx = idx + 1 |
282 |
xfer_bc_tracer(idx) = xfer_array(i,1) |
283 |
ENDDO |
284 |
DO j = 2, Ny |
285 |
idx = idx + 1 |
286 |
xfer_bc_tracer(idx) = xfer_array(Nx,j) |
287 |
ENDDO |
288 |
DO i = (Nx-1), -1, 1 |
289 |
idx = idx + 1 |
290 |
xfer_bc_tracer(idx) = xfer_array(i,Ny) |
291 |
ENDDO |
292 |
DO j = (Ny-1), -1, 2 |
293 |
idx = idx + 1 |
294 |
xfer_bc_tracer(idx) = xfer_array(1,j) |
295 |
ENDDO |
296 |
_BEGIN_MASTER( myThid ) |
297 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
298 |
buffsize = 2*(Nx+Ny)-4 |
299 |
CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION, |
300 |
& local_ice_leader,HsaltBcTag,MPI_COMM_WORLD,mpierr) |
301 |
ENDIF |
302 |
_END_MASTER( myThid ) |
303 |
|
304 |
C Send boundary snow thickness |
305 |
DO bj=1,nSy |
306 |
DO bi=1,nSx |
307 |
DO j=1,sNy |
308 |
DO i=1,sNx |
309 |
local(i,j,bi,bj) = HSNOW(i,j,bi,bj) |
310 |
ENDDO |
311 |
ENDDO |
312 |
ENDDO |
313 |
ENDDO |
314 |
CALL GATHER_2D( xfer_array, local, myThid ) |
315 |
idx = 0 |
316 |
DO i = 1, Nx |
317 |
idx = idx + 1 |
318 |
xfer_bc_tracer(idx) = xfer_array(i,1) |
319 |
ENDDO |
320 |
DO j = 2, Ny |
321 |
idx = idx + 1 |
322 |
xfer_bc_tracer(idx) = xfer_array(Nx,j) |
323 |
ENDDO |
324 |
DO i = (Nx-1), -1, 1 |
325 |
idx = idx + 1 |
326 |
xfer_bc_tracer(idx) = xfer_array(i,Ny) |
327 |
ENDDO |
328 |
DO j = (Ny-1), -1, 2 |
329 |
idx = idx + 1 |
330 |
xfer_bc_tracer(idx) = xfer_array(1,j) |
331 |
ENDDO |
332 |
_BEGIN_MASTER( myThid ) |
333 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
334 |
buffsize = 2*(Nx+Ny)-4 |
335 |
CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION, |
336 |
& local_ice_leader,HsnowBcTag,MPI_COMM_WORLD,mpierr) |
337 |
ENDIF |
338 |
_END_MASTER( myThid ) |
339 |
|
340 |
C Send boundary u ice |
341 |
DO bj=1,nSy |
342 |
DO bi=1,nSx |
343 |
DO j=1,sNy |
344 |
DO i=1,sNx |
345 |
local(i,j,bi,bj) = UICE(i,j,1,bi,bj) |
346 |
ENDDO |
347 |
ENDDO |
348 |
ENDDO |
349 |
ENDDO |
350 |
CALL GATHER_2D( xfer_array, local, myThid ) |
351 |
idx = 0 |
352 |
DO i = 2, Nx |
353 |
idx = idx + 1 |
354 |
xfer_bc_veloc(idx) = xfer_array(i,1) |
355 |
ENDDO |
356 |
DO j = 2, Ny |
357 |
idx = idx + 1 |
358 |
xfer_bc_veloc(idx) = xfer_array(Nx,j) |
359 |
ENDDO |
360 |
DO i = (Nx-1), -1, 2 |
361 |
idx = idx + 1 |
362 |
xfer_bc_veloc(idx) = xfer_array(i,Ny) |
363 |
ENDDO |
364 |
DO j = (Ny-1), -1, 2 |
365 |
idx = idx + 1 |
366 |
xfer_bc_veloc(idx) = xfer_array(2,j) |
367 |
ENDDO |
368 |
_BEGIN_MASTER( myThid ) |
369 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
370 |
buffsize = 2*(Nx+Ny)-6 |
371 |
CALL MPI_SEND(xfer_bc_veloc,buffsize,MPI_DOUBLE_PRECISION, |
372 |
& local_ice_leader,UiceBcTag,MPI_COMM_WORLD,mpierr) |
373 |
ENDIF |
374 |
_END_MASTER( myThid ) |
375 |
|
376 |
C Send boundary v ice |
377 |
DO bj=1,nSy |
378 |
DO bi=1,nSx |
379 |
DO j=1,sNy |
380 |
DO i=1,sNx |
381 |
local(i,j,bi,bj) = VICE(i,j,1,bi,bj) |
382 |
ENDDO |
383 |
ENDDO |
384 |
ENDDO |
385 |
ENDDO |
386 |
CALL GATHER_2D( xfer_array, local, myThid ) |
387 |
idx = 0 |
388 |
DO i = 1, Nx |
389 |
idx = idx + 1 |
390 |
xfer_bc_veloc(idx) = xfer_array(i,2) |
391 |
ENDDO |
392 |
DO j = 3, Ny |
393 |
idx = idx + 1 |
394 |
xfer_bc_veloc(idx) = xfer_array(Nx,j) |
395 |
ENDDO |
396 |
DO i = (Nx-1), -1, 1 |
397 |
idx = idx + 1 |
398 |
xfer_bc_veloc(idx) = xfer_array(i,Ny) |
399 |
ENDDO |
400 |
DO j = (Ny-1), -1, 3 |
401 |
idx = idx + 1 |
402 |
xfer_bc_veloc(idx) = xfer_array(1,j) |
403 |
ENDDO |
404 |
_BEGIN_MASTER( myThid ) |
405 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
406 |
buffsize = 2*(Nx+Ny)-6 |
407 |
CALL MPI_SEND(xfer_bc_veloc,buffsize,MPI_DOUBLE_PRECISION, |
408 |
& local_ice_leader,ViceBcTag,MPI_COMM_WORLD,mpierr) |
409 |
ENDIF |
410 |
_END_MASTER( myThid ) |
411 |
|
412 |
C Send u-wind velocity |
413 |
DO bj=1,nSy |
414 |
DO bi=1,nSx |
415 |
DO j=1,sNy |
416 |
DO i=1,sNx |
417 |
local(i,j,bi,bj) = uwind(i,j,bi,bj) |
418 |
ENDDO |
419 |
ENDDO |
420 |
ENDDO |
421 |
ENDDO |
422 |
CALL GATHER_2D( xfer_array, local, myThid ) |
423 |
_BEGIN_MASTER( myThid ) |
424 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
425 |
#ifdef FIX_FOR_EDGE_WINDS |
426 |
DO j=1,Ny |
427 |
xfer_array(Nx,j)=xfer_array(Nx-1,j) |
428 |
ENDDO |
429 |
#endif |
430 |
buffsize = Nx*Ny |
431 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
432 |
& local_ice_leader,UwindTag,MPI_COMM_WORLD,mpierr) |
433 |
ENDIF |
434 |
_END_MASTER( myThid ) |
435 |
|
436 |
C Send v-wind velocity |
437 |
DO bj=1,nSy |
438 |
DO bi=1,nSx |
439 |
DO j=1,sNy |
440 |
DO i=1,sNx |
441 |
local(i,j,bi,bj) = vwind(i,j,bi,bj) |
442 |
ENDDO |
443 |
ENDDO |
444 |
ENDDO |
445 |
ENDDO |
446 |
CALL GATHER_2D( xfer_array, local, myThid ) |
447 |
_BEGIN_MASTER( myThid ) |
448 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
449 |
#ifdef FIX_FOR_EDGE_WINDS |
450 |
DO i=1,Nx |
451 |
xfer_array(i,Ny)=xfer_array(i,Ny-1) |
452 |
ENDDO |
453 |
#endif |
454 |
buffsize = Nx*Ny |
455 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
456 |
& local_ice_leader,VwindTag,MPI_COMM_WORLD,mpierr) |
457 |
ENDIF |
458 |
_END_MASTER( myThid ) |
459 |
|
460 |
C Send downward longwave radiation |
461 |
DO bj=1,nSy |
462 |
DO bi=1,nSx |
463 |
DO j=1,sNy |
464 |
DO i=1,sNx |
465 |
local(i,j,bi,bj) = lwdown(i,j,bi,bj) |
466 |
ENDDO |
467 |
ENDDO |
468 |
ENDDO |
469 |
ENDDO |
470 |
CALL GATHER_2D( xfer_array, local, myThid ) |
471 |
_BEGIN_MASTER( myThid ) |
472 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
473 |
buffsize = Nx*Ny |
474 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
475 |
& local_ice_leader,LwDownTag,MPI_COMM_WORLD,mpierr) |
476 |
ENDIF |
477 |
_END_MASTER( myThid ) |
478 |
|
479 |
C Send downward shortwave radiation |
480 |
DO bj=1,nSy |
481 |
DO bi=1,nSx |
482 |
DO j=1,sNy |
483 |
DO i=1,sNx |
484 |
local(i,j,bi,bj) = swdown(i,j,bi,bj) |
485 |
ENDDO |
486 |
ENDDO |
487 |
ENDDO |
488 |
ENDDO |
489 |
CALL GATHER_2D( xfer_array, local, myThid ) |
490 |
_BEGIN_MASTER( myThid ) |
491 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
492 |
buffsize = Nx*Ny |
493 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
494 |
& local_ice_leader,SwDownTag,MPI_COMM_WORLD,mpierr) |
495 |
ENDIF |
496 |
_END_MASTER( myThid ) |
497 |
|
498 |
C Send air temperature |
499 |
DO bj=1,nSy |
500 |
DO bi=1,nSx |
501 |
DO j=1,sNy |
502 |
DO i=1,sNx |
503 |
local(i,j,bi,bj) = atemp(i,j,bi,bj) |
504 |
ENDDO |
505 |
ENDDO |
506 |
ENDDO |
507 |
ENDDO |
508 |
CALL GATHER_2D( xfer_array, local, myThid ) |
509 |
_BEGIN_MASTER( myThid ) |
510 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
511 |
buffsize = Nx*Ny |
512 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
513 |
& local_ice_leader,AtempTag,MPI_COMM_WORLD,mpierr) |
514 |
ENDIF |
515 |
_END_MASTER( myThid ) |
516 |
|
517 |
C Send humidity |
518 |
DO bj=1,nSy |
519 |
DO bi=1,nSx |
520 |
DO j=1,sNy |
521 |
DO i=1,sNx |
522 |
local(i,j,bi,bj) = aqh(i,j,bi,bj) |
523 |
ENDDO |
524 |
ENDDO |
525 |
ENDDO |
526 |
ENDDO |
527 |
CALL GATHER_2D( xfer_array, local, myThid ) |
528 |
_BEGIN_MASTER( myThid ) |
529 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
530 |
buffsize = Nx*Ny |
531 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
532 |
& local_ice_leader,AqhTag,MPI_COMM_WORLD,mpierr) |
533 |
ENDIF |
534 |
_END_MASTER( myThid ) |
535 |
|
536 |
C Send precipitation |
537 |
DO bj=1,nSy |
538 |
DO bi=1,nSx |
539 |
DO j=1,sNy |
540 |
DO i=1,sNx |
541 |
local(i,j,bi,bj) = precip(i,j,bi,bj) |
542 |
ENDDO |
543 |
ENDDO |
544 |
ENDDO |
545 |
ENDDO |
546 |
CALL GATHER_2D( xfer_array, local, myThid ) |
547 |
_BEGIN_MASTER( myThid ) |
548 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
549 |
buffsize = Nx*Ny |
550 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
551 |
& local_ice_leader,PrecipTag,MPI_COMM_WORLD,mpierr) |
552 |
ENDIF |
553 |
_END_MASTER( myThid ) |
554 |
|
555 |
C Send ocean surface temperature |
556 |
DO bj=1,nSy |
557 |
DO bi=1,nSx |
558 |
DO j=1,sNy |
559 |
DO i=1,sNx |
560 |
local(i,j,bi,bj) = theta(i,j,1,bi,bj) |
561 |
ENDDO |
562 |
ENDDO |
563 |
ENDDO |
564 |
ENDDO |
565 |
CALL GATHER_2D( xfer_array, local, myThid ) |
566 |
_BEGIN_MASTER( myThid ) |
567 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
568 |
buffsize = Nx*Ny |
569 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
570 |
& local_ice_leader,SstTag,MPI_COMM_WORLD,mpierr) |
571 |
ENDIF |
572 |
_END_MASTER( myThid ) |
573 |
|
574 |
C Send surface u current |
575 |
DO bj=1,nSy |
576 |
DO bi=1,nSx |
577 |
DO j=1,sNy |
578 |
DO i=1,sNx |
579 |
local(i,j,bi,bj) = uVel(i,j,1,bi,bj) |
580 |
ENDDO |
581 |
ENDDO |
582 |
ENDDO |
583 |
ENDDO |
584 |
CALL GATHER_2D( xfer_array, local, myThid ) |
585 |
_BEGIN_MASTER( myThid ) |
586 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
587 |
buffsize = Nx*Ny |
588 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
589 |
& local_ice_leader,UvelTag,MPI_COMM_WORLD,mpierr) |
590 |
ENDIF |
591 |
_END_MASTER( myThid ) |
592 |
|
593 |
C Send surface v current |
594 |
DO bj=1,nSy |
595 |
DO bi=1,nSx |
596 |
DO j=1,sNy |
597 |
DO i=1,sNx |
598 |
local(i,j,bi,bj) = vVel(i,j,1,bi,bj) |
599 |
ENDDO |
600 |
ENDDO |
601 |
ENDDO |
602 |
ENDDO |
603 |
CALL GATHER_2D( xfer_array, local, myThid ) |
604 |
_BEGIN_MASTER( myThid ) |
605 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
606 |
buffsize = Nx*Ny |
607 |
CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
608 |
& local_ice_leader,VvelTag,MPI_COMM_WORLD,mpierr) |
609 |
ENDIF |
610 |
_END_MASTER( myThid ) |
611 |
#ifdef CPL_DEBUG |
612 |
CALL PLOT_FIELD_XYZRL( vVel, 'vVel(k=1)', 1, myIter, myThid ) |
613 |
#endif |
614 |
|
615 |
C Receive ice model time |
616 |
_BEGIN_MASTER( myThid ) |
617 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
618 |
buffsize = 1 |
619 |
CALL MPI_RECV(xfer_scalar,1,MPI_DOUBLE_PRECISION, |
620 |
& local_ice_leader,IceTimeTag,MPI_COMM_WORLD,mpistatus,mpierr) |
621 |
#ifdef CPL_DEBUG |
622 |
print*,'MITgcm receive IceTime', xfer_scalar |
623 |
#endif |
624 |
ENDIF |
625 |
_END_MASTER( myThid ) |
626 |
|
627 |
C Receive ice area Nx*Ny Real*8 |
628 |
_BEGIN_MASTER( myThid ) |
629 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
630 |
buffsize = Nx*Ny |
631 |
CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
632 |
& local_ice_leader,AreaTag,MPI_COMM_WORLD,mpistatus,mpierr) |
633 |
ENDIF |
634 |
_END_MASTER( myThid ) |
635 |
CALL SCATTER_2D( xfer_array, local, myThid ) |
636 |
DO bj=1,nSy |
637 |
DO bi=1,nSx |
638 |
DO j=1,sNy |
639 |
DO i=1,sNx |
640 |
ScatArray(i,j,bi,bj) = local(i,j,bi,bj) |
641 |
ENDDO |
642 |
ENDDO |
643 |
ENDDO |
644 |
ENDDO |
645 |
#ifdef CPL_DEBUG |
646 |
CALL PLOT_FIELD_XYRL( ScatArray, 'ice area', myIter, myThid ) |
647 |
#endif |
648 |
|
649 |
C Receive ice thickness |
650 |
_BEGIN_MASTER( myThid ) |
651 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
652 |
buffsize = Nx*Ny |
653 |
CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
654 |
& local_ice_leader,HeffTag,MPI_COMM_WORLD,mpistatus,mpierr) |
655 |
ENDIF |
656 |
_END_MASTER( myThid ) |
657 |
CALL SCATTER_2D( xfer_array, local, myThid ) |
658 |
DO bj=1,nSy |
659 |
DO bi=1,nSx |
660 |
DO j=1,sNy |
661 |
DO i=1,sNx |
662 |
ScatArray(i,j,bi,bj) = local(i,j,bi,bj) |
663 |
ENDDO |
664 |
ENDDO |
665 |
ENDDO |
666 |
ENDDO |
667 |
#ifdef CPL_DEBUG |
668 |
CALL PLOT_FIELD_XYRL( ScatArray, 'ice thickness', myIter, myThid ) |
669 |
#endif |
670 |
|
671 |
C Receive ice salinity |
672 |
_BEGIN_MASTER( myThid ) |
673 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
674 |
buffsize = Nx*Ny |
675 |
CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
676 |
& local_ice_leader,HsaltTag,MPI_COMM_WORLD,mpistatus,mpierr) |
677 |
ENDIF |
678 |
_END_MASTER( myThid ) |
679 |
CALL SCATTER_2D( xfer_array, local, myThid ) |
680 |
DO bj=1,nSy |
681 |
DO bi=1,nSx |
682 |
DO j=1,sNy |
683 |
DO i=1,sNx |
684 |
ScatArray(i,j,bi,bj) = local(i,j,bi,bj) |
685 |
ENDDO |
686 |
ENDDO |
687 |
ENDDO |
688 |
ENDDO |
689 |
#ifdef CPL_DEBUG |
690 |
CALL PLOT_FIELD_XYRL( ScatArray, 'ice salinity', myIter, myThid ) |
691 |
#endif |
692 |
|
693 |
C Receive snow thickness |
694 |
_BEGIN_MASTER( myThid ) |
695 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
696 |
buffsize = Nx*Ny |
697 |
CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
698 |
& local_ice_leader,HsnowTag,MPI_COMM_WORLD,mpistatus,mpierr) |
699 |
ENDIF |
700 |
_END_MASTER( myThid ) |
701 |
CALL SCATTER_2D( xfer_array, local, myThid ) |
702 |
DO bj=1,nSy |
703 |
DO bi=1,nSx |
704 |
DO j=1,sNy |
705 |
DO i=1,sNx |
706 |
ScatArray(i,j,bi,bj) = local(i,j,bi,bj) |
707 |
ENDDO |
708 |
ENDDO |
709 |
ENDDO |
710 |
ENDDO |
711 |
#ifdef CPL_DEBUG |
712 |
CALL PLOT_FIELD_XYRL( ScatArray, 'ice thickness', myIter, myThid ) |
713 |
#endif |
714 |
|
715 |
C Receive u surface stress |
716 |
_BEGIN_MASTER( myThid ) |
717 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
718 |
buffsize = Nx*Ny |
719 |
CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
720 |
& local_ice_leader,UstressTag,MPI_COMM_WORLD,mpistatus,mpierr) |
721 |
ENDIF |
722 |
_END_MASTER( myThid ) |
723 |
CALL SCATTER_2D( xfer_array, local, myThid ) |
724 |
DO bj=1,nSy |
725 |
DO bi=1,nSx |
726 |
DO j=1,sNy |
727 |
DO i=1,sNx |
728 |
ScatArray(i,j,bi,bj) = local(i,j,bi,bj) |
729 |
ENDDO |
730 |
ENDDO |
731 |
ENDDO |
732 |
ENDDO |
733 |
#ifdef CPL_DEBUG |
734 |
CALL PLOT_FIELD_XYRL( ScatArray, 'u stress', myIter, myThid ) |
735 |
#endif |
736 |
|
737 |
C Receive v surface stress |
738 |
_BEGIN_MASTER( myThid ) |
739 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
740 |
buffsize = Nx*Ny |
741 |
CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
742 |
& local_ice_leader,VstressTag,MPI_COMM_WORLD,mpistatus,mpierr) |
743 |
ENDIF |
744 |
_END_MASTER( myThid ) |
745 |
CALL SCATTER_2D( xfer_array, local, myThid ) |
746 |
DO bj=1,nSy |
747 |
DO bi=1,nSx |
748 |
DO j=1,sNy |
749 |
DO i=1,sNx |
750 |
ScatArray(i,j,bi,bj) = local(i,j,bi,bj) |
751 |
ENDDO |
752 |
ENDDO |
753 |
ENDDO |
754 |
ENDDO |
755 |
#ifdef CPL_DEBUG |
756 |
CALL PLOT_FIELD_XYRL( ScatArray, 'v stress', myIter, myThid ) |
757 |
#endif |
758 |
|
759 |
C Receive residual shortwave |
760 |
_BEGIN_MASTER( myThid ) |
761 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
762 |
buffsize = Nx*Ny |
763 |
CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
764 |
& local_ice_leader,SwResidTag,MPI_COMM_WORLD,mpistatus,mpierr) |
765 |
ENDIF |
766 |
_END_MASTER( myThid ) |
767 |
CALL SCATTER_2D( xfer_array, local, myThid ) |
768 |
DO bj=1,nSy |
769 |
DO bi=1,nSx |
770 |
DO j=1,sNy |
771 |
DO i=1,sNx |
772 |
ScatArray(i,j,bi,bj) = local(i,j,bi,bj) |
773 |
ENDDO |
774 |
ENDDO |
775 |
ENDDO |
776 |
ENDDO |
777 |
#ifdef CPL_DEBUG |
778 |
CALL PLOT_FIELD_XYRL( ScatArray, 'shortwave', myIter, myThid ) |
779 |
#endif |
780 |
|
781 |
C Receive heat flux |
782 |
_BEGIN_MASTER( myThid ) |
783 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
784 |
buffsize = Nx*Ny |
785 |
CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
786 |
& local_ice_leader,HeatFluxTag,MPI_COMM_WORLD,mpistatus,mpierr) |
787 |
ENDIF |
788 |
_END_MASTER( myThid ) |
789 |
CALL SCATTER_2D( xfer_array, local, myThid ) |
790 |
DO bj=1,nSy |
791 |
DO bi=1,nSx |
792 |
DO j=1,sNy |
793 |
DO i=1,sNx |
794 |
ScatArray(i,j,bi,bj) = local(i,j,bi,bj) |
795 |
ENDDO |
796 |
ENDDO |
797 |
ENDDO |
798 |
ENDDO |
799 |
#ifdef CPL_DEBUG |
800 |
CALL PLOT_FIELD_XYRL( ScatArray, 'heat flux', myIter, myThid ) |
801 |
#endif |
802 |
|
803 |
C Receive freshwater flux |
804 |
_BEGIN_MASTER( myThid ) |
805 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
806 |
buffsize = Nx*Ny |
807 |
CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
808 |
& local_ice_leader,WaterFluxTag,MPI_COMM_WORLD,mpistatus,mpierr) |
809 |
ENDIF |
810 |
_END_MASTER( myThid ) |
811 |
CALL SCATTER_2D( xfer_array, local, myThid ) |
812 |
DO bj=1,nSy |
813 |
DO bi=1,nSx |
814 |
DO j=1,sNy |
815 |
DO i=1,sNx |
816 |
ScatArray(i,j,bi,bj) = local(i,j,bi,bj) |
817 |
ENDDO |
818 |
ENDDO |
819 |
ENDDO |
820 |
ENDDO |
821 |
#ifdef CPL_DEBUG |
822 |
CALL PLOT_FIELD_XYRL( ScatArray, 'freshwater', myIter, myThid ) |
823 |
#endif |
824 |
|
825 |
C Receive salt flux |
826 |
_BEGIN_MASTER( myThid ) |
827 |
IF ( myworldid .EQ. local_ocean_leader ) THEN |
828 |
buffsize = Nx*Ny |
829 |
CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, |
830 |
& local_ice_leader,SaltFluxTag,MPI_COMM_WORLD,mpistatus,mpierr) |
831 |
ENDIF |
832 |
_END_MASTER( myThid ) |
833 |
CALL SCATTER_2D( xfer_array, local, myThid ) |
834 |
DO bj=1,nSy |
835 |
DO bi=1,nSx |
836 |
DO j=1,sNy |
837 |
DO i=1,sNx |
838 |
ScatArray(i,j,bi,bj) = local(i,j,bi,bj) |
839 |
ENDDO |
840 |
ENDDO |
841 |
ENDDO |
842 |
ENDDO |
843 |
#ifdef CPL_DEBUG |
844 |
CALL PLOT_FIELD_XYRL( ScatArray, 'salt flux', myIter, myThid ) |
845 |
#endif |
846 |
|
847 |
#endif /* ALLOW_CPL_MPMICE */ |
848 |
|
849 |
RETURN |
850 |
END |