/[MITgcm]/MITgcm_contrib/MPMice/beaufort/code/cpl_mpmice.F
ViewVC logotype

Contents of /MITgcm_contrib/MPMice/beaufort/code/cpl_mpmice.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.16 - (show annotations) (download)
Sun Mar 18 01:25:21 2012 UTC (13 years, 4 months ago) by dimitri
Branch: MAIN
Changes since 1.15: +1 -1 lines
moving itername inside CPL_DEBUG

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

  ViewVC Help
Powered by ViewVC 1.1.22