/[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.20 - (show annotations) (download)
Wed Sep 28 15:11:47 2016 UTC (8 years, 10 months ago) by dimitri
Branch: MAIN
CVS Tags: HEAD
Changes since 1.19: +1 -0 lines
Fixing bug that made ice melt too much during summer months.
I was double-counting shortwave radiation.

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 longitude East of center of grid cell
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) = xC(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( xC, 'xC', 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,xCTag,MPI_COMM_WORLD,mpierr)
120 ENDIF
121 _END_MASTER( myThid )
122 # endif /* CPL_COUPLED */
123
124 C Send latitude North of center of grid cell
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) = yC(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( yC, 'yC', 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,yCTag,MPI_COMM_WORLD,mpierr)
144 ENDIF
145 _END_MASTER( myThid )
146 # endif /* CPL_COUPLED */
147
148 C Send longitude East of SouthWest corner
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) = xG(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( xG, 'xG', 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,xGTag,MPI_COMM_WORLD,mpierr)
168 ENDIF
169 _END_MASTER( myThid )
170 # endif /* CPL_COUPLED */
171
172 C Send latitude North of SouthWest corner
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) = yG(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( yG, 'yG', 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,yGTag,MPI_COMM_WORLD,mpierr)
192 ENDIF
193 _END_MASTER( myThid )
194 # endif /* CPL_COUPLED */
195
196 C Send distance in m between SouthWest and SouthEast corner
197 DO bj=1,nSy
198 DO bi=1,nSx
199 DO j=1,sNy
200 DO i=1,sNx
201 local(i,j,bi,bj) = dxG(i,j,bi,bj)
202 ENDDO
203 ENDDO
204 ENDDO
205 ENDDO
206 CALL GATHER_2D( xfer_array, local, myThid )
207 # ifdef CPL_DEBUG
208 CALL PLOT_FIELD_XYRL( dxG, 'dxG', myIter, myThid )
209 # endif /* CPL_DEBUG */
210 # ifdef CPL_COUPLED
211 _BEGIN_MASTER( myThid )
212 IF ( myworldid .EQ. local_ocean_leader ) THEN
213 buffsize = Nx*Ny
214 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
215 & local_ice_leader,dxGTag,MPI_COMM_WORLD,mpierr)
216 ENDIF
217 _END_MASTER( myThid )
218 # endif /* CPL_COUPLED */
219
220 C Send distance in m between SouthWest and NorthEast corner
221 DO bj=1,nSy
222 DO bi=1,nSx
223 DO j=1,sNy
224 DO i=1,sNx
225 local(i,j,bi,bj) = dyG(i,j,bi,bj)
226 ENDDO
227 ENDDO
228 ENDDO
229 ENDDO
230 CALL GATHER_2D( xfer_array, local, myThid )
231 # ifdef CPL_DEBUG
232 CALL PLOT_FIELD_XYRL( dyG, 'dyG', myIter, myThid )
233 # endif /* CPL_DEBUG */
234 # ifdef CPL_COUPLED
235 _BEGIN_MASTER( myThid )
236 IF ( myworldid .EQ. local_ocean_leader ) THEN
237 buffsize = Nx*Ny
238 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
239 & local_ice_leader,dyGTag,MPI_COMM_WORLD,mpierr)
240 ENDIF
241 _END_MASTER( myThid )
242 # endif /* CPL_COUPLED */
243
244 C Send cosine(alpha) relative to geographic direction at grid cell center
245 DO bj=1,nSy
246 DO bi=1,nSx
247 DO j=1,sNy
248 DO i=1,sNx
249 local(i,j,bi,bj) = angleCosC(i,j,bi,bj)
250 ENDDO
251 ENDDO
252 ENDDO
253 ENDDO
254 CALL GATHER_2D( xfer_array, local, myThid )
255 # ifdef CPL_DEBUG
256 CALL PLOT_FIELD_XYRL( angleCosC, 'aCS', myIter, myThid )
257 # endif /* CPL_DEBUG */
258 # ifdef CPL_COUPLED
259 _BEGIN_MASTER( myThid )
260 IF ( myworldid .EQ. local_ocean_leader ) THEN
261 buffsize = Nx*Ny
262 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
263 & local_ice_leader,aCStag,MPI_COMM_WORLD,mpierr)
264 ENDIF
265 _END_MASTER( myThid )
266 # endif /* CPL_COUPLED */
267
268 C Send sine(alpha) relative to geographic direction at grid cell center
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) = angleSinC(i,j,bi,bj)
274 ENDDO
275 ENDDO
276 ENDDO
277 ENDDO
278 CALL GATHER_2D( xfer_array, local, myThid )
279 # ifdef CPL_DEBUG
280 CALL PLOT_FIELD_XYRL( angleSinC, 'aSN', myIter, myThid )
281 # endif /* CPL_DEBUG */
282 # ifdef CPL_COUPLED
283 _BEGIN_MASTER( myThid )
284 IF ( myworldid .EQ. local_ocean_leader ) THEN
285 buffsize = Nx*Ny
286 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
287 & local_ice_leader,aSNtag,MPI_COMM_WORLD,mpierr)
288 ENDIF
289 _END_MASTER( myThid )
290 # endif /* CPL_COUPLED */
291
292 C Send landmask of center of grid cell, 0 is land, >0 is ocean
293 DO bj=1,nSy
294 DO bi=1,nSx
295 DO j=1,sNy
296 DO i=1,sNx
297 local(i,j,bi,bj) = hFacC(i,j,1,bi,bj)
298 ENDDO
299 ENDDO
300 ENDDO
301 ENDDO
302 CALL GATHER_2D( xfer_array, local, myThid )
303 # ifdef CPL_DEBUG
304 CALL PLOT_FIELD_XYRL( hFacC, 'hFacC', myIter, myThid )
305 # endif /* CPL_DEBUG */
306 # ifdef CPL_COUPLED
307 _BEGIN_MASTER( myThid )
308 IF ( myworldid .EQ. local_ocean_leader ) THEN
309 buffsize = Nx*Ny
310 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
311 & local_ice_leader,hFacCtag,MPI_COMM_WORLD,mpierr)
312 ENDIF
313 _END_MASTER( myThid )
314 # endif /* CPL_COUPLED */
315
316 C Send ice area
317 DO bj=1,nSy
318 DO bi=1,nSx
319 DO j=1,sNy
320 DO i=1,sNx
321 local(i,j,bi,bj) = AREA(i,j,bi,bj)
322 ENDDO
323 ENDDO
324 ENDDO
325 ENDDO
326 CALL GATHER_2D( xfer_array, local, myThid )
327 # ifdef CPL_DEBUG
328 CALL PLOT_FIELD_XYRL( AREA, 'AREA', myIter, myThid )
329 # endif /* CPL_DEBUG */
330 # ifdef CPL_COUPLED
331 _BEGIN_MASTER( myThid )
332 IF ( myworldid .EQ. local_ocean_leader ) THEN
333 buffsize = Nx*Ny
334 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
335 & local_ice_leader,AreaTag,MPI_COMM_WORLD,mpierr)
336 ENDIF
337 _END_MASTER( myThid )
338 # endif /* CPL_COUPLED */
339
340 C Send ice thickness
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) = HEFF(i,j,bi,bj)
346 ENDDO
347 ENDDO
348 ENDDO
349 ENDDO
350 CALL GATHER_2D( xfer_array, local, myThid )
351 # ifdef CPL_DEBUG
352 CALL PLOT_FIELD_XYRL( HEFF, 'HEFF', myIter, myThid )
353 # endif /* CPL_DEBUG */
354 # ifdef CPL_COUPLED
355 _BEGIN_MASTER( myThid )
356 IF ( myworldid .EQ. local_ocean_leader ) THEN
357 buffsize = Nx*Ny
358 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
359 & local_ice_leader,HeffTag,MPI_COMM_WORLD,mpierr)
360 ENDIF
361 _END_MASTER( myThid )
362 # endif /* CPL_COUPLED */
363
364 C Send ice salinity
365 DO bj=1,nSy
366 DO bi=1,nSx
367 DO j=1,sNy
368 DO i=1,sNx
369 local(i,j,bi,bj) = HSALT(i,j,bi,bj)
370 ENDDO
371 ENDDO
372 ENDDO
373 ENDDO
374 CALL GATHER_2D( xfer_array, local, myThid )
375 # ifdef CPL_DEBUG
376 CALL PLOT_FIELD_XYRL( HSALT, 'HSALT', myIter, myThid )
377 # endif /* CPL_DEBUG */
378 # ifdef CPL_COUPLED
379 _BEGIN_MASTER( myThid )
380 IF ( myworldid .EQ. local_ocean_leader ) THEN
381 buffsize = Nx*Ny
382 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
383 & local_ice_leader,HsaltTag,MPI_COMM_WORLD,mpierr)
384 ENDIF
385 _END_MASTER( myThid )
386 # endif /* CPL_COUPLED */
387
388 C Send snow thickness
389 DO bj=1,nSy
390 DO bi=1,nSx
391 DO j=1,sNy
392 DO i=1,sNx
393 local(i,j,bi,bj) = HSNOW(i,j,bi,bj)
394 ENDDO
395 ENDDO
396 ENDDO
397 ENDDO
398 CALL GATHER_2D( xfer_array, local, myThid )
399 # ifdef CPL_DEBUG
400 CALL PLOT_FIELD_XYRL( HSNOW, 'HSNOW', myIter, myThid )
401 # endif /* CPL_DEBUG */
402 # ifdef CPL_COUPLED
403 _BEGIN_MASTER( myThid )
404 IF ( myworldid .EQ. local_ocean_leader ) THEN
405 buffsize = Nx*Ny
406 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
407 & local_ice_leader,HsnowTag,MPI_COMM_WORLD,mpierr)
408 ENDIF
409 _END_MASTER( myThid )
410 # endif /* CPL_COUPLED */
411
412 ENDIF ! ( myTime .EQ. startTime )
413
414 C-- Apply ice open boundary conditions
415 #ifdef ALLOW_OBCS
416 IF ( useOBCS ) THEN
417 CALL OBCS_APPLY_SEAICE( myThid )
418 CALL OBCS_APPLY_UVICE( uice, vice, myThid )
419 ENDIF
420 #endif /* ALLOW_OBCS */
421
422 C Send ocean model time
423 _BEGIN_MASTER( myThid )
424 xfer_scalar = myTime
425 buffsize = 1
426 # ifdef CPL_DEBUG
427 print*,'MITgcm send OceanTime', xfer_scalar
428 # endif /* CPL_DEBUG */
429 # ifdef CPL_COUPLED
430 IF ( myworldid .EQ. local_ocean_leader ) THEN
431 CALL MPI_SEND(xfer_scalar,buffsize,MPI_DOUBLE_PRECISION,
432 & local_ice_leader,OceanTimeTag,MPI_COMM_WORLD,mpierr)
433 ENDIF
434 # endif /* CPL_COUPLED */
435 _END_MASTER( myThid )
436
437 C Send boundary ice area
438 DO bj=1,nSy
439 DO bi=1,nSx
440 DO j=1,sNy
441 DO i=1,sNx
442 local(i,j,bi,bj) = AREA(i,j,bi,bj)
443 ENDDO
444 ENDDO
445 ENDDO
446 ENDDO
447 CALL GATHER_2D( xfer_array, local, myThid )
448 idx = 0
449 DO i = 1, Nx
450 idx = idx + 1
451 xfer_bc_tracer(idx) = xfer_array(i,1)
452 ENDDO
453 DO j = 2, Ny
454 idx = idx + 1
455 xfer_bc_tracer(idx) = xfer_array(Nx,j)
456 ENDDO
457 DO i = (Nx-1), 1, -1
458 idx = idx + 1
459 xfer_bc_tracer(idx) = xfer_array(i,Ny)
460 ENDDO
461 DO j = (Ny-1), 2, -1
462 idx = idx + 1
463 xfer_bc_tracer(idx) = xfer_array(1,j)
464 ENDDO
465 buffsize = 2*(Nx+Ny)-4
466 # ifdef CPL_DEBUG
467 CALL PLOT_FIELD_XYRL( AREA, 'AREA obcs', myIter, myThid )
468 CALL WRITE_GLVEC_RS ( 'AREAobcs.', itername,
469 & xfer_bc_tracer, buffsize, myIter, myThid )
470 # endif /* CPL_DEBUG */
471 # ifdef CPL_COUPLED
472 _BEGIN_MASTER( myThid )
473 IF ( myworldid .EQ. local_ocean_leader ) THEN
474 cdb print*,'MITgcm is about to send AreaBcTag',buffsize
475 CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION,
476 & local_ice_leader,AreaBcTag,MPI_COMM_WORLD,mpierr)
477 cdb print*,'MITgcm has sent AreaBcTag',buffsize
478 ENDIF
479 _END_MASTER( myThid )
480 # endif /* CPL_COUPLED */
481
482 C Send boundary ice thickness
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) = HEFF(i,j,bi,bj)
488 ENDDO
489 ENDDO
490 ENDDO
491 ENDDO
492 CALL GATHER_2D( xfer_array, local, myThid )
493 idx = 0
494 DO i = 1, Nx
495 idx = idx + 1
496 xfer_bc_tracer(idx) = xfer_array(i,1)
497 ENDDO
498 DO j = 2, Ny
499 idx = idx + 1
500 xfer_bc_tracer(idx) = xfer_array(Nx,j)
501 ENDDO
502 DO i = (Nx-1), 1, -1
503 idx = idx + 1
504 xfer_bc_tracer(idx) = xfer_array(i,Ny)
505 ENDDO
506 DO j = (Ny-1), 2, -1
507 idx = idx + 1
508 xfer_bc_tracer(idx) = xfer_array(1,j)
509 ENDDO
510 buffsize = 2*(Nx+Ny)-4
511 # ifdef CPL_DEBUG
512 CALL PLOT_FIELD_XYRL( HEFF, 'HEFF obcs', myIter, myThid )
513 CALL WRITE_GLVEC_RS ( 'HEFFobcs.', itername,
514 & xfer_bc_tracer, buffsize, myIter, myThid )
515 # endif /* CPL_DEBUG */
516 # ifdef CPL_COUPLED
517 _BEGIN_MASTER( myThid )
518 IF ( myworldid .EQ. local_ocean_leader ) THEN
519 CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION,
520 & local_ice_leader,HeffBcTag,MPI_COMM_WORLD,mpierr)
521 ENDIF
522 _END_MASTER( myThid )
523 # endif /* CPL_COUPLED */
524
525 C Send boundary ice salinity
526 DO bj=1,nSy
527 DO bi=1,nSx
528 DO j=1,sNy
529 DO i=1,sNx
530 local(i,j,bi,bj) = HSALT(i,j,bi,bj)
531 ENDDO
532 ENDDO
533 ENDDO
534 ENDDO
535 CALL GATHER_2D( xfer_array, local, myThid )
536 idx = 0
537 DO i = 1, Nx
538 idx = idx + 1
539 xfer_bc_tracer(idx) = xfer_array(i,1)
540 ENDDO
541 DO j = 2, Ny
542 idx = idx + 1
543 xfer_bc_tracer(idx) = xfer_array(Nx,j)
544 ENDDO
545 DO i = (Nx-1), 1, -1
546 idx = idx + 1
547 xfer_bc_tracer(idx) = xfer_array(i,Ny)
548 ENDDO
549 DO j = (Ny-1), 2, -1
550 idx = idx + 1
551 xfer_bc_tracer(idx) = xfer_array(1,j)
552 ENDDO
553 buffsize = 2*(Nx+Ny)-4
554 # ifdef CPL_DEBUG
555 CALL PLOT_FIELD_XYRL( HSALT, 'HSALT obcs', myIter, myThid )
556 CALL WRITE_GLVEC_RS ( 'HSALTobcs.', itername,
557 & xfer_bc_tracer, buffsize, myIter, myThid )
558 # endif /* CPL_DEBUG */
559 # ifdef CPL_COUPLED
560 _BEGIN_MASTER( myThid )
561 IF ( myworldid .EQ. local_ocean_leader ) THEN
562 CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION,
563 & local_ice_leader,HsaltBcTag,MPI_COMM_WORLD,mpierr)
564 ENDIF
565 _END_MASTER( myThid )
566 # endif /* CPL_COUPLED */
567
568 C Send boundary snow thickness
569 DO bj=1,nSy
570 DO bi=1,nSx
571 DO j=1,sNy
572 DO i=1,sNx
573 local(i,j,bi,bj) = HSNOW(i,j,bi,bj)
574 ENDDO
575 ENDDO
576 ENDDO
577 ENDDO
578 CALL GATHER_2D( xfer_array, local, myThid )
579 idx = 0
580 DO i = 1, Nx
581 idx = idx + 1
582 xfer_bc_tracer(idx) = xfer_array(i,1)
583 ENDDO
584 DO j = 2, Ny
585 idx = idx + 1
586 xfer_bc_tracer(idx) = xfer_array(Nx,j)
587 ENDDO
588 DO i = (Nx-1), 1, -1
589 idx = idx + 1
590 xfer_bc_tracer(idx) = xfer_array(i,Ny)
591 ENDDO
592 DO j = (Ny-1), 2, -1
593 idx = idx + 1
594 xfer_bc_tracer(idx) = xfer_array(1,j)
595 ENDDO
596 buffsize = 2*(Nx+Ny)-4
597 # ifdef CPL_DEBUG
598 CALL PLOT_FIELD_XYRL( HSNOW, 'HSNOW obcs', myIter, myThid )
599 CALL WRITE_GLVEC_RS ( 'HSNOWobcs.', itername,
600 & xfer_bc_tracer, buffsize, myIter, myThid )
601 # endif /* CPL_DEBUG */
602 # ifdef CPL_COUPLED
603 _BEGIN_MASTER( myThid )
604 IF ( myworldid .EQ. local_ocean_leader ) THEN
605 CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION,
606 & local_ice_leader,HsnowBcTag,MPI_COMM_WORLD,mpierr)
607 ENDIF
608 _END_MASTER( myThid )
609 # endif /* CPL_COUPLED */
610
611 C Send boundary u ice
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) = UICE(i,j,bi,bj)
617 ENDDO
618 ENDDO
619 ENDDO
620 ENDDO
621 CALL GATHER_2D( xfer_array, local, myThid )
622 idx = 0
623 DO i = 2, Nx
624 idx = idx + 1
625 xfer_bc_veloc(idx) = xfer_array(i,1)
626 ENDDO
627 DO j = 2, Ny
628 idx = idx + 1
629 xfer_bc_veloc(idx) = xfer_array(Nx,j)
630 ENDDO
631 DO i = (Nx-1), 2, -1
632 idx = idx + 1
633 xfer_bc_veloc(idx) = xfer_array(i,Ny)
634 ENDDO
635 DO j = (Ny-1), 2, -1
636 idx = idx + 1
637 xfer_bc_veloc(idx) = xfer_array(2,j)
638 ENDDO
639 buffsize = 2*(Nx+Ny)-6
640 # ifdef CPL_DEBUG
641 CALL PLOT_FIELD_XYRL( UICE, 'UICE obcs', myIter, myThid )
642 CALL WRITE_GLVEC_RS ( 'UICEobcs.', itername,
643 & xfer_bc_veloc, buffsize, myIter, myThid )
644 # endif /* CPL_DEBUG */
645 # ifdef CPL_COUPLED
646 _BEGIN_MASTER( myThid )
647 IF ( myworldid .EQ. local_ocean_leader ) THEN
648 CALL MPI_SEND(xfer_bc_veloc,buffsize,MPI_DOUBLE_PRECISION,
649 & local_ice_leader,UiceBcTag,MPI_COMM_WORLD,mpierr)
650 ENDIF
651 _END_MASTER( myThid )
652 # endif /* CPL_COUPLED */
653
654 C Send boundary v ice
655 DO bj=1,nSy
656 DO bi=1,nSx
657 DO j=1,sNy
658 DO i=1,sNx
659 local(i,j,bi,bj) = VICE(i,j,bi,bj)
660 ENDDO
661 ENDDO
662 ENDDO
663 ENDDO
664 CALL GATHER_2D( xfer_array, local, myThid )
665 idx = 0
666 DO i = 1, Nx
667 idx = idx + 1
668 xfer_bc_veloc(idx) = xfer_array(i,2)
669 ENDDO
670 DO j = 3, Ny
671 idx = idx + 1
672 xfer_bc_veloc(idx) = xfer_array(Nx,j)
673 ENDDO
674 DO i = (Nx-1), 1, -1
675 idx = idx + 1
676 xfer_bc_veloc(idx) = xfer_array(i,Ny)
677 ENDDO
678 DO j = (Ny-1), 3, -1
679 idx = idx + 1
680 xfer_bc_veloc(idx) = xfer_array(1,j)
681 ENDDO
682 buffsize = 2*(Nx+Ny)-6
683 # ifdef CPL_DEBUG
684 CALL PLOT_FIELD_XYRL( VICE, 'VICE obcs', myIter, myThid )
685 CALL WRITE_GLVEC_RS ( 'VICEobcs.', itername,
686 & xfer_bc_veloc, buffsize, myIter, myThid )
687 # endif /* CPL_DEBUG */
688 # ifdef CPL_COUPLED
689 _BEGIN_MASTER( myThid )
690 IF ( myworldid .EQ. local_ocean_leader ) THEN
691 CALL MPI_SEND(xfer_bc_veloc,buffsize,MPI_DOUBLE_PRECISION,
692 & local_ice_leader,ViceBcTag,MPI_COMM_WORLD,mpierr)
693 ENDIF
694 _END_MASTER( myThid )
695 # endif /* CPL_COUPLED */
696
697 C Send u-wind velocity
698 DO bj=1,nSy
699 DO bi=1,nSx
700 DO j=1,sNy
701 DO i=1,sNx
702 local(i,j,bi,bj) = uwind(i,j,bi,bj)
703 ENDDO
704 ENDDO
705 ENDDO
706 ENDDO
707 CALL GATHER_2D( xfer_array, local, myThid )
708 # ifdef CPL_DEBUG
709 CALL PLOT_FIELD_XYRL( UWIND, 'UWIND', myIter, myThid )
710 # endif /* CPL_DEBUG */
711 # ifdef CPL_COUPLED
712 _BEGIN_MASTER( myThid )
713 IF ( myworldid .EQ. local_ocean_leader ) THEN
714 # ifdef FIX_FOR_EDGE_WINDS
715 DO j=1,Ny
716 xfer_array(Nx,j)=xfer_array(Nx-1,j)
717 ENDDO
718 # endif /* FIX_FOR_EDGE_WINDS */
719 buffsize = Nx*Ny
720 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
721 & local_ice_leader,UwindTag,MPI_COMM_WORLD,mpierr)
722 ENDIF
723 _END_MASTER( myThid )
724 # endif /* CPL_COUPLED */
725
726 C Send v-wind velocity
727 DO bj=1,nSy
728 DO bi=1,nSx
729 DO j=1,sNy
730 DO i=1,sNx
731 local(i,j,bi,bj) = vwind(i,j,bi,bj)
732 ENDDO
733 ENDDO
734 ENDDO
735 ENDDO
736 CALL GATHER_2D( xfer_array, local, myThid )
737 # ifdef CPL_DEBUG
738 CALL PLOT_FIELD_XYRL( VWIND, 'VWIND', myIter, myThid )
739 # endif /* CPL_DEBUG */
740 # ifdef CPL_COUPLED
741 _BEGIN_MASTER( myThid )
742 IF ( myworldid .EQ. local_ocean_leader ) THEN
743 # ifdef FIX_FOR_EDGE_WINDS
744 DO i=1,Nx
745 xfer_array(i,Ny)=xfer_array(i,Ny-1)
746 ENDDO
747 # endif /* FIX_FOR_EDGE_WINDS */
748 buffsize = Nx*Ny
749 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
750 & local_ice_leader,VwindTag,MPI_COMM_WORLD,mpierr)
751 ENDIF
752 _END_MASTER( myThid )
753 # endif /* CPL_COUPLED */
754
755 C Send downward longwave radiation
756 DO bj=1,nSy
757 DO bi=1,nSx
758 DO j=1,sNy
759 DO i=1,sNx
760 local(i,j,bi,bj) = lwdown(i,j,bi,bj)
761 ENDDO
762 ENDDO
763 ENDDO
764 ENDDO
765 CALL GATHER_2D( xfer_array, local, myThid )
766 # ifdef CPL_DEBUG
767 CALL PLOT_FIELD_XYRL( LWDOWN, 'LWDOWN', myIter, myThid )
768 # endif /* CPL_DEBUG */
769 # ifdef CPL_COUPLED
770 _BEGIN_MASTER( myThid )
771 IF ( myworldid .EQ. local_ocean_leader ) THEN
772 buffsize = Nx*Ny
773 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
774 & local_ice_leader,LwDownTag,MPI_COMM_WORLD,mpierr)
775 ENDIF
776 _END_MASTER( myThid )
777 # endif /* CPL_COUPLED */
778
779 C Send downward shortwave radiation
780 DO bj=1,nSy
781 DO bi=1,nSx
782 DO j=1,sNy
783 DO i=1,sNx
784 local(i,j,bi,bj) = swdown(i,j,bi,bj)
785 ENDDO
786 ENDDO
787 ENDDO
788 ENDDO
789 CALL GATHER_2D( xfer_array, local, myThid )
790 # ifdef CPL_DEBUG
791 CALL PLOT_FIELD_XYRL( SWDOWN, 'SWDOWN', myIter, myThid )
792 # endif /* CPL_DEBUG */
793 # ifdef CPL_COUPLED
794 _BEGIN_MASTER( myThid )
795 IF ( myworldid .EQ. local_ocean_leader ) THEN
796 buffsize = Nx*Ny
797 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
798 & local_ice_leader,SwDownTag,MPI_COMM_WORLD,mpierr)
799 ENDIF
800 _END_MASTER( myThid )
801 # endif /* CPL_COUPLED */
802
803 C Send air temperature
804 DO bj=1,nSy
805 DO bi=1,nSx
806 DO j=1,sNy
807 DO i=1,sNx
808 local(i,j,bi,bj) = atemp(i,j,bi,bj)
809 ENDDO
810 ENDDO
811 ENDDO
812 ENDDO
813 CALL GATHER_2D( xfer_array, local, myThid )
814 # ifdef CPL_DEBUG
815 CALL PLOT_FIELD_XYRL( ATEMP, 'ATEMP', myIter, myThid )
816 # endif /* CPL_DEBUG */
817 # ifdef CPL_COUPLED
818 _BEGIN_MASTER( myThid )
819 IF ( myworldid .EQ. local_ocean_leader ) THEN
820 buffsize = Nx*Ny
821 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
822 & local_ice_leader,AtempTag,MPI_COMM_WORLD,mpierr)
823 ENDIF
824 _END_MASTER( myThid )
825 # endif /* CPL_COUPLED */
826
827 C Send humidity
828 DO bj=1,nSy
829 DO bi=1,nSx
830 DO j=1,sNy
831 DO i=1,sNx
832 local(i,j,bi,bj) = aqh(i,j,bi,bj)
833 ENDDO
834 ENDDO
835 ENDDO
836 ENDDO
837 CALL GATHER_2D( xfer_array, local, myThid )
838 # ifdef CPL_DEBUG
839 CALL PLOT_FIELD_XYRL( AQH, 'AQH', myIter, myThid )
840 # endif /* CPL_DEBUG */
841 # ifdef CPL_COUPLED
842 _BEGIN_MASTER( myThid )
843 IF ( myworldid .EQ. local_ocean_leader ) THEN
844 buffsize = Nx*Ny
845 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
846 & local_ice_leader,AqhTag,MPI_COMM_WORLD,mpierr)
847 ENDIF
848 _END_MASTER( myThid )
849 # endif /* CPL_COUPLED */
850
851 C Send precipitation
852 DO bj=1,nSy
853 DO bi=1,nSx
854 DO j=1,sNy
855 DO i=1,sNx
856 local(i,j,bi,bj) = precip(i,j,bi,bj)
857 ENDDO
858 ENDDO
859 ENDDO
860 ENDDO
861 CALL GATHER_2D( xfer_array, local, myThid )
862 # ifdef CPL_DEBUG
863 CALL PLOT_FIELD_XYRL( PRECIP, 'PRECIP', myIter, myThid )
864 # endif /* CPL_DEBUG */
865 # ifdef CPL_COUPLED
866 _BEGIN_MASTER( myThid )
867 IF ( myworldid .EQ. local_ocean_leader ) THEN
868 buffsize = Nx*Ny
869 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
870 & local_ice_leader,PrecipTag,MPI_COMM_WORLD,mpierr)
871 ENDIF
872 _END_MASTER( myThid )
873 # endif /* CPL_COUPLED */
874
875 C Send ocean surface temperature
876 DO bj=1,nSy
877 DO bi=1,nSx
878 DO j=1,sNy
879 DO i=1,sNx
880 local(i,j,bi,bj) = theta(i,j,1,bi,bj)
881 ENDDO
882 ENDDO
883 ENDDO
884 ENDDO
885 CALL GATHER_2D( xfer_array, local, myThid )
886 # ifdef CPL_DEBUG
887 CALL PLOT_FIELD_XYZRL( THETA, 'SST', 1, myIter, myThid )
888 # endif /* CPL_DEBUG */
889 # ifdef CPL_COUPLED
890 _BEGIN_MASTER( myThid )
891 IF ( myworldid .EQ. local_ocean_leader ) THEN
892 buffsize = Nx*Ny
893 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
894 & local_ice_leader,SstTag,MPI_COMM_WORLD,mpierr)
895 ENDIF
896 _END_MASTER( myThid )
897 # endif /* CPL_COUPLED */
898
899 C Send ocean surface salinity
900 DO bj=1,nSy
901 DO bi=1,nSx
902 DO j=1,sNy
903 DO i=1,sNx
904 local(i,j,bi,bj) = salt(i,j,1,bi,bj)
905 ENDDO
906 ENDDO
907 ENDDO
908 ENDDO
909 CALL GATHER_2D( xfer_array, local, myThid )
910 # ifdef CPL_DEBUG
911 CALL PLOT_FIELD_XYZRL( SALT, 'SSS', 1, myIter, myThid )
912 # endif /* CPL_DEBUG */
913 # ifdef CPL_COUPLED
914 _BEGIN_MASTER( myThid )
915 IF ( myworldid .EQ. local_ocean_leader ) THEN
916 buffsize = Nx*Ny
917 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
918 & local_ice_leader,SssTag,MPI_COMM_WORLD,mpierr)
919 ENDIF
920 _END_MASTER( myThid )
921 # endif /* CPL_COUPLED */
922
923 C Send surface u current
924 DO bj=1,nSy
925 DO bi=1,nSx
926 DO j=1,sNy
927 DO i=1,sNx
928 local(i,j,bi,bj) = uVel(i,j,1,bi,bj)
929 ENDDO
930 ENDDO
931 ENDDO
932 ENDDO
933 CALL GATHER_2D( xfer_array, local, myThid )
934 # ifdef CPL_DEBUG
935 CALL PLOT_FIELD_XYZRL( uVel, 'uVel(k=1)', 1, myIter, myThid )
936 # endif /* CPL_DEBUG */
937 # ifdef CPL_COUPLED
938 _BEGIN_MASTER( myThid )
939 IF ( myworldid .EQ. local_ocean_leader ) THEN
940 buffsize = Nx*Ny
941 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
942 & local_ice_leader,UvelTag,MPI_COMM_WORLD,mpierr)
943 ENDIF
944 _END_MASTER( myThid )
945 # endif /* CPL_COUPLED */
946
947 C Send surface v current
948 DO bj=1,nSy
949 DO bi=1,nSx
950 DO j=1,sNy
951 DO i=1,sNx
952 local(i,j,bi,bj) = vVel(i,j,1,bi,bj)
953 ENDDO
954 ENDDO
955 ENDDO
956 ENDDO
957 CALL GATHER_2D( xfer_array, local, myThid )
958 # ifdef CPL_DEBUG
959 CALL PLOT_FIELD_XYZRL( vVel, 'vVel(k=1)', 1, myIter, myThid )
960 # endif /* CPL_DEBUG */
961 # ifdef CPL_COUPLED
962 _BEGIN_MASTER( myThid )
963 IF ( myworldid .EQ. local_ocean_leader ) THEN
964 buffsize = Nx*Ny
965 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
966 & local_ice_leader,VvelTag,MPI_COMM_WORLD,mpierr)
967 ENDIF
968 _END_MASTER( myThid )
969 # endif /* CPL_COUPLED */
970
971 C Receive ice model time
972 _BEGIN_MASTER( myThid )
973 # ifdef CPL_DEBUG
974 print*,'MITgcm receive IceTime'
975 # endif /* CPL_DEBUG */
976 # ifdef CPL_COUPLED
977 IF ( myworldid .EQ. local_ocean_leader ) THEN
978 buffsize = 1
979 CALL MPI_RECV(xfer_scalar,1,MPI_DOUBLE_PRECISION,
980 & local_ice_leader,IceTimeTag,MPI_COMM_WORLD,mpistatus,mpierr)
981 ENDIF
982 # endif /* CPL_COUPLED */
983 _END_MASTER( myThid )
984
985 C Receive ice area Nx*Ny Real*8
986 # ifdef CPL_COUPLED
987 _BEGIN_MASTER( myThid )
988 IF ( myworldid .EQ. local_ocean_leader ) THEN
989 buffsize = Nx*Ny
990 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
991 & local_ice_leader,AreaTag,MPI_COMM_WORLD,mpistatus,mpierr)
992 ENDIF
993 _END_MASTER( myThid )
994 CALL SCATTER_2D( xfer_array, local, myThid )
995 DO bj=1,nSy
996 DO bi=1,nSx
997 DO j=1,sNy
998 DO i=1,sNx
999 AREA(i,j,bi,bj) = local(i,j,bi,bj)
1000 ENDDO
1001 ENDDO
1002 ENDDO
1003 ENDDO
1004 # endif /* CPL_COUPLED */
1005 # ifdef CPL_DEBUG
1006 CALL PLOT_FIELD_XYRL( AREA, 'ice area', myIter, myThid )
1007 # endif /* CPL_DEBUG */
1008
1009 C Receive ice thickness
1010 # ifdef CPL_COUPLED
1011 _BEGIN_MASTER( myThid )
1012 IF ( myworldid .EQ. local_ocean_leader ) THEN
1013 buffsize = Nx*Ny
1014 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
1015 & local_ice_leader,HeffTag,MPI_COMM_WORLD,mpistatus,mpierr)
1016 ENDIF
1017 _END_MASTER( myThid )
1018 CALL SCATTER_2D( xfer_array, local, myThid )
1019 DO bj=1,nSy
1020 DO bi=1,nSx
1021 DO j=1,sNy
1022 DO i=1,sNx
1023 HEFF(i,j,bi,bj) = local(i,j,bi,bj)
1024 ENDDO
1025 ENDDO
1026 ENDDO
1027 ENDDO
1028 # endif /* CPL_COUPLED */
1029 # ifdef CPL_DEBUG
1030 CALL PLOT_FIELD_XYRL( HEFF, 'ice thickness', myIter, myThid )
1031 # endif /* CPL_DEBUG */
1032
1033 C Receive ice salinity
1034 # ifdef CPL_COUPLED
1035 _BEGIN_MASTER( myThid )
1036 IF ( myworldid .EQ. local_ocean_leader ) THEN
1037 buffsize = Nx*Ny
1038 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
1039 & local_ice_leader,HsaltTag,MPI_COMM_WORLD,mpistatus,mpierr)
1040 ENDIF
1041 _END_MASTER( myThid )
1042 CALL SCATTER_2D( xfer_array, local, myThid )
1043 DO bj=1,nSy
1044 DO bi=1,nSx
1045 DO j=1,sNy
1046 DO i=1,sNx
1047 HSALT(i,j,bi,bj) = local(i,j,bi,bj)
1048 ENDDO
1049 ENDDO
1050 ENDDO
1051 ENDDO
1052 # endif /* CPL_COUPLED */
1053 # ifdef CPL_DEBUG
1054 CALL PLOT_FIELD_XYRL( HSALT, 'ice salinity', myIter, myThid )
1055 # endif /* CPL_DEBUG */
1056
1057 C Receive snow thickness
1058 # ifdef CPL_COUPLED
1059 _BEGIN_MASTER( myThid )
1060 IF ( myworldid .EQ. local_ocean_leader ) THEN
1061 buffsize = Nx*Ny
1062 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
1063 & local_ice_leader,HsnowTag,MPI_COMM_WORLD,mpistatus,mpierr)
1064 ENDIF
1065 _END_MASTER( myThid )
1066 CALL SCATTER_2D( xfer_array, local, myThid )
1067 DO bj=1,nSy
1068 DO bi=1,nSx
1069 DO j=1,sNy
1070 DO i=1,sNx
1071 HSNOW(i,j,bi,bj) = local(i,j,bi,bj)
1072 ENDDO
1073 ENDDO
1074 ENDDO
1075 ENDDO
1076 # endif /* CPL_COUPLED */
1077 # ifdef CPL_DEBUG
1078 CALL PLOT_FIELD_XYRL( HSNOW, 'snow thickness', myIter, myThid )
1079 # endif /* CPL_DEBUG */
1080
1081 C Receive u ice velocity
1082 # ifdef CPL_COUPLED
1083 _BEGIN_MASTER( myThid )
1084 IF ( myworldid .EQ. local_ocean_leader ) THEN
1085 buffsize = Nx*Ny
1086 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
1087 & local_ice_leader,UiceTag,MPI_COMM_WORLD,mpistatus,mpierr)
1088 ENDIF
1089 _END_MASTER( myThid )
1090 CALL SCATTER_2D( xfer_array, local, myThid )
1091 DO bj=1,nSy
1092 DO bi=1,nSx
1093 DO j=1,sNy
1094 DO i=1,sNx
1095 UICE(i,j,bi,bj) = local(i,j,bi,bj)
1096 ENDDO
1097 ENDDO
1098 ENDDO
1099 ENDDO
1100 # ifdef CPL_DEBUG
1101 CALL PLOT_FIELD_XYRL( local, 'uice', myIter, myThid )
1102 # endif /* CPL_DEBUG */
1103 # endif /* CPL_COUPLED */
1104 # ifdef CPL_DEBUG
1105 CALL PLOT_FIELD_XYRL( UICE, 'uice', myIter, myThid )
1106 # endif /* CPL_DEBUG */
1107
1108 C Receive v ice velocity
1109 # ifdef CPL_COUPLED
1110 _BEGIN_MASTER( myThid )
1111 IF ( myworldid .EQ. local_ocean_leader ) THEN
1112 buffsize = Nx*Ny
1113 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
1114 & local_ice_leader,ViceTag,MPI_COMM_WORLD,mpistatus,mpierr)
1115 ENDIF
1116 _END_MASTER( myThid )
1117 CALL SCATTER_2D( xfer_array, local, myThid )
1118 DO bj=1,nSy
1119 DO bi=1,nSx
1120 DO j=1,sNy
1121 DO i=1,sNx
1122 VICE(i,j,bi,bj) = local(i,j,bi,bj)
1123 ENDDO
1124 ENDDO
1125 ENDDO
1126 ENDDO
1127 # ifdef CPL_DEBUG
1128 CALL PLOT_FIELD_XYRL( local, 'vice', myIter, myThid )
1129 # endif /* CPL_DEBUG */
1130 # endif /* CPL_COUPLED */
1131 # ifdef CPL_DEBUG
1132 CALL PLOT_FIELD_XYRL( VICE, 'vice', myIter, myThid )
1133 # endif /* CPL_DEBUG */
1134
1135 C Receive u surface stress
1136 # ifdef CPL_COUPLED
1137 _BEGIN_MASTER( myThid )
1138 IF ( myworldid .EQ. local_ocean_leader ) THEN
1139 buffsize = Nx*Ny
1140 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
1141 & local_ice_leader,UstressTag,MPI_COMM_WORLD,mpistatus,mpierr)
1142 ENDIF
1143 _END_MASTER( myThid )
1144 CALL SCATTER_2D( xfer_array, local, myThid )
1145 DO bj=1,nSy
1146 DO bi=1,nSx
1147 DO j=1,sNy
1148 DO i=1,sNx
1149 fu(i,j,bi,bj) = AREA(i,j,bi,bj) * local(i,j,bi,bj) +
1150 & (1.-AREA(i,j,bi,bj)) * fu (i,j,bi,bj)
1151 ENDDO
1152 ENDDO
1153 ENDDO
1154 ENDDO
1155 # ifdef CPL_DEBUG
1156 CALL PLOT_FIELD_XYRL( local, 'mpm u stress', myIter, myThid )
1157 # endif /* CPL_DEBUG */
1158 # endif /* CPL_COUPLED */
1159 # ifdef CPL_DEBUG
1160 CALL PLOT_FIELD_XYRL( fu, 'u stress', myIter, myThid )
1161 # endif /* CPL_DEBUG */
1162
1163 C Receive v surface stress
1164 # ifdef CPL_COUPLED
1165 _BEGIN_MASTER( myThid )
1166 IF ( myworldid .EQ. local_ocean_leader ) THEN
1167 buffsize = Nx*Ny
1168 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
1169 & local_ice_leader,VstressTag,MPI_COMM_WORLD,mpistatus,mpierr)
1170 ENDIF
1171 _END_MASTER( myThid )
1172 CALL SCATTER_2D( xfer_array, local, myThid )
1173 DO bj=1,nSy
1174 DO bi=1,nSx
1175 DO j=1,sNy
1176 DO i=1,sNx
1177 fv(i,j,bi,bj) = AREA(i,j,bi,bj) * local(i,j,bi,bj) +
1178 & (1.-AREA(i,j,bi,bj)) * fv (i,j,bi,bj)
1179 ENDDO
1180 ENDDO
1181 ENDDO
1182 ENDDO
1183 # ifdef CPL_DEBUG
1184 CALL PLOT_FIELD_XYRL( local, 'mpm v stress', myIter, myThid )
1185 # endif /* CPL_DEBUG */
1186 # endif /* CPL_COUPLED */
1187 # ifdef CPL_DEBUG
1188 CALL PLOT_FIELD_XYRL( fv, 'v stress', myIter, myThid )
1189 # endif /* CPL_DEBUG */
1190
1191 C Receive residual shortwave
1192 # ifdef CPL_COUPLED
1193 _BEGIN_MASTER( myThid )
1194 IF ( myworldid .EQ. local_ocean_leader ) THEN
1195 buffsize = Nx*Ny
1196 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
1197 & local_ice_leader,SwResidTag,MPI_COMM_WORLD,mpistatus,mpierr)
1198 ENDIF
1199 _END_MASTER( myThid )
1200 CALL SCATTER_2D( xfer_array, local, myThid )
1201 DO bj=1,nSy
1202 DO bi=1,nSx
1203 DO j=1,sNy
1204 DO i=1,sNx
1205 Qnet(i,j,bi,bj) = Qnet(i,j,bi,bj) - Qsw(i,j,bi,bj)
1206 Qsw(i,j,bi,bj) = -AREA(i,j,bi,bj) * local(i,j,bi,bj) +
1207 & (1.-AREA(i,j,bi,bj)) * Qsw(i,j,bi,bj)
1208 ENDDO
1209 ENDDO
1210 ENDDO
1211 ENDDO
1212 # ifdef CPL_DEBUG
1213 CALL PLOT_FIELD_XYRL( local, 'mpm shortwave', myIter, myThid )
1214 # endif /* CPL_DEBUG */
1215 # endif /* CPL_COUPLED */
1216 # ifdef CPL_DEBUG
1217 CALL PLOT_FIELD_XYRL( Qsw, 'shortwave', myIter, myThid )
1218 # endif /* CPL_DEBUG */
1219
1220 C Receive heat flux
1221 # ifdef CPL_COUPLED
1222 _BEGIN_MASTER( myThid )
1223 IF ( myworldid .EQ. local_ocean_leader ) THEN
1224 buffsize = Nx*Ny
1225 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
1226 & local_ice_leader,HeatFluxTag,MPI_COMM_WORLD,mpistatus,mpierr)
1227 ENDIF
1228 _END_MASTER( myThid )
1229 CALL SCATTER_2D( xfer_array, local, myThid )
1230 DO bj=1,nSy
1231 DO bi=1,nSx
1232 DO j=1,sNy
1233 DO i=1,sNx
1234 Qnet(i,j,bi,bj) = Qsw(i,j,bi,bj) -
1235 & AREA(i,j,bi,bj) * local(i,j,bi,bj) +
1236 & (1.-AREA(i,j,bi,bj)) * Qnet(i,j,bi,bj)
1237 ENDDO
1238 ENDDO
1239 ENDDO
1240 ENDDO
1241 # ifdef CPL_DEBUG
1242 CALL PLOT_FIELD_XYRL( local, 'mpm heat flux', myIter, myThid )
1243 # endif /* CPL_DEBUG */
1244 # endif /* CPL_COUPLED */
1245 # ifdef CPL_DEBUG
1246 CALL PLOT_FIELD_XYRL( Qnet, 'heat flux', myIter, myThid )
1247 # endif /* CPL_DEBUG */
1248
1249 C Receive freshwater flux
1250 # ifdef CPL_COUPLED
1251 _BEGIN_MASTER( myThid )
1252 IF ( myworldid .EQ. local_ocean_leader ) THEN
1253 buffsize = Nx*Ny
1254 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
1255 & local_ice_leader,WaterFluxTag,MPI_COMM_WORLD,mpistatus,mpierr)
1256 ENDIF
1257 _END_MASTER( myThid )
1258 CALL SCATTER_2D( xfer_array, local, myThid )
1259 DO bj=1,nSy
1260 DO bi=1,nSx
1261 DO j=1,sNy
1262 DO i=1,sNx
1263 EmPmR(i,j,bi,bj) = - AREA(i,j,bi,bj) * local(i,j,bi,bj) +
1264 & ( 1. - AREA(i,j,bi,bj)) * EmPmR(i,j,bi,bj)
1265 ENDDO
1266 ENDDO
1267 ENDDO
1268 ENDDO
1269 # ifdef CPL_DEBUG
1270 CALL PLOT_FIELD_XYRL( local, 'mpm freshwater', myIter, myThid )
1271 # endif /* CPL_DEBUG */
1272 # endif /* CPL_COUPLED */
1273 # ifdef CPL_DEBUG
1274 CALL PLOT_FIELD_XYRL( EmPmR, 'freshwater', myIter, myThid )
1275 # endif /* CPL_DEBUG */
1276
1277 C Receive salt flux
1278 # ifdef CPL_COUPLED
1279 _BEGIN_MASTER( myThid )
1280 IF ( myworldid .EQ. local_ocean_leader ) THEN
1281 buffsize = Nx*Ny
1282 CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
1283 & local_ice_leader,SaltFluxTag,MPI_COMM_WORLD,mpistatus,mpierr)
1284 ENDIF
1285 _END_MASTER( myThid )
1286 CALL SCATTER_2D( xfer_array, local, myThid )
1287 DO bj=1,nSy
1288 DO bi=1,nSx
1289 DO j=1,sNy
1290 DO i=1,sNx
1291 saltFlux(i,j,bi,bj) = - AREA(i,j,bi,bj) * local(i,j,bi,bj)
1292 ENDDO
1293 ENDDO
1294 ENDDO
1295 ENDDO
1296 # ifdef CPL_DEBUG
1297 CALL PLOT_FIELD_XYRL( local, 'mpm salt flux', myIter, myThid )
1298 # endif /* CPL_DEBUG */
1299 # endif /* CPL_COUPLED */
1300 # ifdef CPL_DEBUG
1301 CALL PLOT_FIELD_XYRL( saltFlux, 'salt flux', myIter, myThid )
1302 # endif /* CPL_DEBUG */
1303
1304 #endif /* ALLOW_CPL_MPMICE */
1305
1306 RETURN
1307 END

  ViewVC Help
Powered by ViewVC 1.1.22