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

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

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


Revision 1.2 - (hide annotations) (download)
Mon Jun 22 19:40:11 2009 UTC (16 years, 1 month ago) by dimitri
Branch: MAIN
Changes since 1.1: +24 -12 lines
Initializing ScatArray to avoid crashing some compilers.
Problem reported by Oksana Guba.

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

  ViewVC Help
Powered by ViewVC 1.1.22