/[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.14 - (show annotations) (download)
Thu Mar 15 20:02:56 2012 UTC (13 years, 4 months ago) by dimitri
Branch: MAIN
Changes since 1.13: +23 -6 lines
adding obcs diagnostics

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

  ViewVC Help
Powered by ViewVC 1.1.22