/[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.9 - (show annotations) (download)
Fri Feb 3 19:13:16 2012 UTC (13 years, 6 months ago) by dimitri
Branch: MAIN
Changes since 1.8: +1 -0 lines
first running two-way coupled configuration on lozenge

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

  ViewVC Help
Powered by ViewVC 1.1.22