/[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.8 - (hide annotations) (download)
Fri Feb 3 18:12:59 2012 UTC (13 years, 5 months ago) by dimitri
Branch: MAIN
Changes since 1.7: +6 -7 lines
more fixes

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

  ViewVC Help
Powered by ViewVC 1.1.22