#define FIX_FOR_EDGE_WINDS #include "PACKAGES_CONFIG.h" #include "CPP_OPTIONS.h" CBOP C !ROUTINE: CPL_MPMICE C !INTERFACE: SUBROUTINE CPL_MPMICE( myTime, myIter, myThid ) C !DESCRIPTION: \bv C *================================================================== C | SUBROUTINE cpl_mpmice C | o Couple MITgcm ocean model with MPMice sea ice model C *================================================================== C \ev C !USES: IMPLICIT NONE C == Global variables == #include "SIZE.h" #include "EEPARAMS.h" #include "PARAMS.h" #include "DYNVARS.h" #include "GRID.h" #include "FFIELDS.h" #include "SEAICE_OPTIONS.h" #include "SEAICE_SIZE.h" #include "SEAICE.h" #ifdef ALLOW_EXF # include "EXF_OPTIONS.h" # include "EXF_FIELDS.h" #endif LOGICAL DIFFERENT_MULTIPLE EXTERNAL DIFFERENT_MULTIPLE C !LOCAL VARIABLES: C mytime - time counter for this thread (seconds) C myiter - iteration counter for this thread C mythid - thread number for this instance of the routine. _RL mytime INTEGER myiter, mythid CEOP #ifdef ALLOW_CPL_MPMICE # include "EESUPPORT.h" # include "CPL_MPMICE.h" COMMON /CPL_MPI_ID/ & myworldid, local_ocean_leader, local_ice_leader integer myworldid, local_ocean_leader, local_ice_leader # ifdef ALLOW_USE_MPI integer mpistatus(MPI_STATUS_SIZE), mpierr # endif /* ALLOW_USE_MPI */ integer xfer_gridsize(2) integer i, j, bi, bj, buffsize, idx Real*8 xfer_scalar Real*8 xfer_array(Nx,Ny) Real*8 xfer_bc_tracer(2*(Nx+Ny)-4) Real*8 xfer_bc_veloc(2*(Nx+Ny)-6) _RL local(1:sNx,1:sNy,nSx,nSy) # ifdef CPL_DEBUG character*(10) itername write(itername,'(i10.10)') myIter # endif /* CPL_DEBUG */ IF( myTime .EQ. startTime ) THEN C Send deltatimestep _BEGIN_MASTER( myThid ) xfer_scalar = deltat buffsize = 1 # ifdef CPL_DEBUG print*,'MITgcm send TimeInterval', xfer_scalar # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED IF ( myworldid .EQ. local_ocean_leader ) THEN CALL MPI_SEND(xfer_scalar,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,TimeIntervalTag,MPI_COMM_WORLD,mpierr) ENDIF # endif /* CPL_COUPLED */ _END_MASTER( myThid ) C Send grid dimensions (Nx,Ny) _BEGIN_MASTER( myThid ) xfer_gridsize(1)=Nx xfer_gridsize(2)=Ny buffsize = 2 # ifdef CPL_DEBUG print*,'MITgcm send OceanGridsize', xfer_gridsize # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED IF ( myworldid .EQ. local_ocean_leader ) THEN CALL MPI_SEND(xfer_gridsize,buffsize,MPI_INTEGER, & local_ice_leader,OceanGridsizeTag,MPI_COMM_WORLD,mpierr) ENDIF # endif /* CPL_COUPLED */ _END_MASTER( myThid ) C Send longitude East of center of grid cell DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = xC(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( xC, 'xC', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,xCTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send latitude North of center of grid cell DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = yC(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( yC, 'yC', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,yCTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send longitude East of SouthWest corner DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = xG(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( xG, 'xG', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,xGTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send latitude North of SouthWest corner DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = yG(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( yG, 'yG', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,yGTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send distance in m between SouthWest and SouthEast corner DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = dxG(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( dxG, 'dxG', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,dxGTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send distance in m between SouthWest and NorthEast corner DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = dyG(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( dyG, 'dyG', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,dyGTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send cosine(alpha) relative to geographic direction at grid cell center DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = angleCosC(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( angleCosC, 'aCS', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,aCStag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send sine(alpha) relative to geographic direction at grid cell center DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = angleSinC(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( angleSinC, 'aSN', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,aSNtag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send landmask of center of grid cell, 0 is land, >0 is ocean DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = hFacC(i,j,1,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( hFacC, 'hFacC', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,hFacCtag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send ice area DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = AREA(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( AREA, 'AREA', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,AreaTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send ice thickness DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = HEFF(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( HEFF, 'HEFF', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,HeffTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send ice salinity DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = HSALT(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( HSALT, 'HSALT', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,HsaltTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send snow thickness DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = HSNOW(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( HSNOW, 'HSNOW', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,HsnowTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ ENDIF ! ( myTime .EQ. startTime ) C-- Apply ice open boundary conditions #ifdef ALLOW_OBCS IF ( useOBCS ) THEN CALL OBCS_APPLY_SEAICE( myThid ) CALL OBCS_APPLY_UVICE( uice, vice, myThid ) ENDIF #endif /* ALLOW_OBCS */ C Send ocean model time _BEGIN_MASTER( myThid ) xfer_scalar = myTime buffsize = 1 # ifdef CPL_DEBUG print*,'MITgcm send OceanTime', xfer_scalar # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED IF ( myworldid .EQ. local_ocean_leader ) THEN CALL MPI_SEND(xfer_scalar,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,OceanTimeTag,MPI_COMM_WORLD,mpierr) ENDIF # endif /* CPL_COUPLED */ _END_MASTER( myThid ) C Send boundary ice area DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = AREA(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) idx = 0 DO i = 1, Nx idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(i,1) ENDDO DO j = 2, Ny idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(Nx,j) ENDDO DO i = (Nx-1), 1, -1 idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(i,Ny) ENDDO DO j = (Ny-1), 2, -1 idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(1,j) ENDDO buffsize = 2*(Nx+Ny)-4 # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( AREA, 'AREA obcs', myIter, myThid ) CALL WRITE_GLVEC_RS ( 'AREAobcs.', itername, & xfer_bc_tracer, buffsize, myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN cdb print*,'MITgcm is about to send AreaBcTag',buffsize CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,AreaBcTag,MPI_COMM_WORLD,mpierr) cdb print*,'MITgcm has sent AreaBcTag',buffsize ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send boundary ice thickness DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = HEFF(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) idx = 0 DO i = 1, Nx idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(i,1) ENDDO DO j = 2, Ny idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(Nx,j) ENDDO DO i = (Nx-1), 1, -1 idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(i,Ny) ENDDO DO j = (Ny-1), 2, -1 idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(1,j) ENDDO buffsize = 2*(Nx+Ny)-4 # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( HEFF, 'HEFF obcs', myIter, myThid ) CALL WRITE_GLVEC_RS ( 'HEFFobcs.', itername, & xfer_bc_tracer, buffsize, myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,HeffBcTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send boundary ice salinity DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = HSALT(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) idx = 0 DO i = 1, Nx idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(i,1) ENDDO DO j = 2, Ny idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(Nx,j) ENDDO DO i = (Nx-1), 1, -1 idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(i,Ny) ENDDO DO j = (Ny-1), 2, -1 idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(1,j) ENDDO buffsize = 2*(Nx+Ny)-4 # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( HSALT, 'HSALT obcs', myIter, myThid ) CALL WRITE_GLVEC_RS ( 'HSALTobcs.', itername, & xfer_bc_tracer, buffsize, myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,HsaltBcTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send boundary snow thickness DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = HSNOW(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) idx = 0 DO i = 1, Nx idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(i,1) ENDDO DO j = 2, Ny idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(Nx,j) ENDDO DO i = (Nx-1), 1, -1 idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(i,Ny) ENDDO DO j = (Ny-1), 2, -1 idx = idx + 1 xfer_bc_tracer(idx) = xfer_array(1,j) ENDDO buffsize = 2*(Nx+Ny)-4 # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( HSNOW, 'HSNOW obcs', myIter, myThid ) CALL WRITE_GLVEC_RS ( 'HSNOWobcs.', itername, & xfer_bc_tracer, buffsize, myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN CALL MPI_SEND(xfer_bc_tracer,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,HsnowBcTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send boundary u ice DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = UICE(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) idx = 0 DO i = 2, Nx idx = idx + 1 xfer_bc_veloc(idx) = xfer_array(i,1) ENDDO DO j = 2, Ny idx = idx + 1 xfer_bc_veloc(idx) = xfer_array(Nx,j) ENDDO DO i = (Nx-1), 2, -1 idx = idx + 1 xfer_bc_veloc(idx) = xfer_array(i,Ny) ENDDO DO j = (Ny-1), 2, -1 idx = idx + 1 xfer_bc_veloc(idx) = xfer_array(2,j) ENDDO buffsize = 2*(Nx+Ny)-6 # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( UICE, 'UICE obcs', myIter, myThid ) CALL WRITE_GLVEC_RS ( 'UICEobcs.', itername, & xfer_bc_veloc, buffsize, myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN CALL MPI_SEND(xfer_bc_veloc,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,UiceBcTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send boundary v ice DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = VICE(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) idx = 0 DO i = 1, Nx idx = idx + 1 xfer_bc_veloc(idx) = xfer_array(i,2) ENDDO DO j = 3, Ny idx = idx + 1 xfer_bc_veloc(idx) = xfer_array(Nx,j) ENDDO DO i = (Nx-1), 1, -1 idx = idx + 1 xfer_bc_veloc(idx) = xfer_array(i,Ny) ENDDO DO j = (Ny-1), 3, -1 idx = idx + 1 xfer_bc_veloc(idx) = xfer_array(1,j) ENDDO buffsize = 2*(Nx+Ny)-6 # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( VICE, 'VICE obcs', myIter, myThid ) CALL WRITE_GLVEC_RS ( 'VICEobcs.', itername, & xfer_bc_veloc, buffsize, myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN CALL MPI_SEND(xfer_bc_veloc,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,ViceBcTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send u-wind velocity DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = uwind(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( UWIND, 'UWIND', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN # ifdef FIX_FOR_EDGE_WINDS DO j=1,Ny xfer_array(Nx,j)=xfer_array(Nx-1,j) ENDDO # endif /* FIX_FOR_EDGE_WINDS */ buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,UwindTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send v-wind velocity DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = vwind(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( VWIND, 'VWIND', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN # ifdef FIX_FOR_EDGE_WINDS DO i=1,Nx xfer_array(i,Ny)=xfer_array(i,Ny-1) ENDDO # endif /* FIX_FOR_EDGE_WINDS */ buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,VwindTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send downward longwave radiation DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = lwdown(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( LWDOWN, 'LWDOWN', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,LwDownTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send downward shortwave radiation DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = swdown(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( SWDOWN, 'SWDOWN', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,SwDownTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send air temperature DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = atemp(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( ATEMP, 'ATEMP', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,AtempTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send humidity DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = aqh(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( AQH, 'AQH', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,AqhTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send precipitation DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = precip(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( PRECIP, 'PRECIP', myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,PrecipTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send ocean surface temperature DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = theta(i,j,1,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYZRL( THETA, 'SST', 1, myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,SstTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send ocean surface salinity DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = salt(i,j,1,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYZRL( SALT, 'SSS', 1, myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,SssTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send surface u current DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = uVel(i,j,1,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYZRL( uVel, 'uVel(k=1)', 1, myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,UvelTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Send surface v current DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx local(i,j,bi,bj) = vVel(i,j,1,bi,bj) ENDDO ENDDO ENDDO ENDDO CALL GATHER_2D( xfer_array, local, myThid ) # ifdef CPL_DEBUG CALL PLOT_FIELD_XYZRL( vVel, 'vVel(k=1)', 1, myIter, myThid ) # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,VvelTag,MPI_COMM_WORLD,mpierr) ENDIF _END_MASTER( myThid ) # endif /* CPL_COUPLED */ C Receive ice model time _BEGIN_MASTER( myThid ) # ifdef CPL_DEBUG print*,'MITgcm receive IceTime' # endif /* CPL_DEBUG */ # ifdef CPL_COUPLED IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = 1 CALL MPI_RECV(xfer_scalar,1,MPI_DOUBLE_PRECISION, & local_ice_leader,IceTimeTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF # endif /* CPL_COUPLED */ _END_MASTER( myThid ) C Receive ice area Nx*Ny Real*8 # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,AreaTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx AREA(i,j,bi,bj) = local(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( AREA, 'ice area', myIter, myThid ) # endif /* CPL_DEBUG */ C Receive ice thickness # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,HeffTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx HEFF(i,j,bi,bj) = local(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( HEFF, 'ice thickness', myIter, myThid ) # endif /* CPL_DEBUG */ C Receive ice salinity # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,HsaltTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx HSALT(i,j,bi,bj) = local(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( HSALT, 'ice salinity', myIter, myThid ) # endif /* CPL_DEBUG */ C Receive snow thickness # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,HsnowTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx HSNOW(i,j,bi,bj) = local(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( HSNOW, 'snow thickness', myIter, myThid ) # endif /* CPL_DEBUG */ C Receive u ice velocity # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,UiceTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx UICE(i,j,bi,bj) = local(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( local, 'uice', myIter, myThid ) # endif /* CPL_DEBUG */ # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( UICE, 'uice', myIter, myThid ) # endif /* CPL_DEBUG */ C Receive v ice velocity # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,ViceTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx VICE(i,j,bi,bj) = local(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( local, 'vice', myIter, myThid ) # endif /* CPL_DEBUG */ # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( VICE, 'vice', myIter, myThid ) # endif /* CPL_DEBUG */ C Receive u surface stress # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,UstressTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx fu(i,j,bi,bj) = AREA(i,j,bi,bj) * local(i,j,bi,bj) + & (1.-AREA(i,j,bi,bj)) * fu (i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( local, 'mpm u stress', myIter, myThid ) # endif /* CPL_DEBUG */ # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( fu, 'u stress', myIter, myThid ) # endif /* CPL_DEBUG */ C Receive v surface stress # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,VstressTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx fv(i,j,bi,bj) = AREA(i,j,bi,bj) * local(i,j,bi,bj) + & (1.-AREA(i,j,bi,bj)) * fv (i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( local, 'mpm v stress', myIter, myThid ) # endif /* CPL_DEBUG */ # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( fv, 'v stress', myIter, myThid ) # endif /* CPL_DEBUG */ C Receive residual shortwave # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,SwResidTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx Qnet(i,j,bi,bj) = Qnet(i,j,bi,bj) - Qsw(i,j,bi,bj) Qsw(i,j,bi,bj) = -AREA(i,j,bi,bj) * local(i,j,bi,bj) + & (1.-AREA(i,j,bi,bj)) * Qsw(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( local, 'mpm shortwave', myIter, myThid ) # endif /* CPL_DEBUG */ # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( Qsw, 'shortwave', myIter, myThid ) # endif /* CPL_DEBUG */ C Receive heat flux # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,HeatFluxTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx Qnet(i,j,bi,bj) = Qsw(i,j,bi,bj) - & AREA(i,j,bi,bj) * local(i,j,bi,bj) + & (1.-AREA(i,j,bi,bj)) * Qnet(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( local, 'mpm heat flux', myIter, myThid ) # endif /* CPL_DEBUG */ # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( Qnet, 'heat flux', myIter, myThid ) # endif /* CPL_DEBUG */ C Receive freshwater flux # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,WaterFluxTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx EmPmR(i,j,bi,bj) = - AREA(i,j,bi,bj) * local(i,j,bi,bj) + & ( 1. - AREA(i,j,bi,bj)) * EmPmR(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( local, 'mpm freshwater', myIter, myThid ) # endif /* CPL_DEBUG */ # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( EmPmR, 'freshwater', myIter, myThid ) # endif /* CPL_DEBUG */ C Receive salt flux # ifdef CPL_COUPLED _BEGIN_MASTER( myThid ) IF ( myworldid .EQ. local_ocean_leader ) THEN buffsize = Nx*Ny CALL MPI_RECV(xfer_array,buffsize,MPI_DOUBLE_PRECISION, & local_ice_leader,SaltFluxTag,MPI_COMM_WORLD,mpistatus,mpierr) ENDIF _END_MASTER( myThid ) CALL SCATTER_2D( xfer_array, local, myThid ) DO bj=1,nSy DO bi=1,nSx DO j=1,sNy DO i=1,sNx saltFlux(i,j,bi,bj) = - AREA(i,j,bi,bj) * local(i,j,bi,bj) ENDDO ENDDO ENDDO ENDDO # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( local, 'mpm salt flux', myIter, myThid ) # endif /* CPL_DEBUG */ # endif /* CPL_COUPLED */ # ifdef CPL_DEBUG CALL PLOT_FIELD_XYRL( saltFlux, 'salt flux', myIter, myThid ) # endif /* CPL_DEBUG */ #endif /* ALLOW_CPL_MPMICE */ RETURN END