Changeset 22546


Ignore:
Timestamp:
03/16/18 12:09:04 (7 years ago)
Author:
dmenemen
Message:

MITgcm sends XC/YC to ISSM.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • issm/trunk-jpl/test/MITgcm/code/cpl_issm.F

    r22532 r22546  
    4848     &     mpiMyWid, toissmcomm
    4949      INTEGER mpiMyWid, toissmcomm
    50       integer mpistatus(MPI_STATUS_SIZE), mpiRC
     50      integer mpistatus(MPI_STATUS_SIZE)
    5151      real*8 CouplingTime, IceModelTime
     52      integer i, j, bi, bj, buffsize
     53      Real*8  xfer_array(Nx,Ny)
     54      _RL     local(1:sNx,1:sNy,nSx,nSy)
    5255
    53       print*,'choubi 1'
    54      
    5556      IF( myTime .EQ. startTime ) THEN
    56       print*,'choubi 2'
    5757
    5858C     Send deltatimestep
    59        _BEGIN_MASTER( myThid )
    60 
    61        print*,'choubi 3'
    62 
     59         _BEGIN_MASTER( myThid )
    6360         call MPI_Recv(CouplingTime,1,MPI_DOUBLE,0,10001000,toissmcomm,
    64      &        mpistatus,mpiRC)
     61     &        mpistatus)
    6562         print*, 'Ocean received CouplingTime: ', CouplingTime
    6663         call MPI_Recv(IceModelTime,1,MPI_DOUBLE,0,10001001,toissmcomm,
    67      &        mpistatus,mpiRC)
     64     &        mpistatus)
    6865         print*, 'Ocean received IceModelTime: ', IceModelTime
    6966         call MPI_Send(myTime,1,MPI_DOUBLE,0,10001002,toissmcomm,
    70      &        mpistatus,mpiRC)
     67     &        mpistatus)
    7168         call MPI_Send(Nx,1,MPI_INT,0,10001003,toissmcomm,
    72      &        mpistatus,mpiRC)
     69     &        mpistatus)
    7370         call MPI_Send(Ny,1,MPI_INT,0,10001004,toissmcomm,
    74      &        mpistatus,mpiRC)
    75      
    76 cdb      call MPI_Bcast(dummy1,1,MPI_INT,0,MPI_COMM_MODEL,mpiRC)
    77 cdb      print*, 'Ocean : dummy received:',dummy1
     71     &        mpistatus)
     72         _END_MASTER( myThid )
    7873
    79       _END_MASTER( myThid )
     74C     Send longitude East of center of cell
     75         DO bj=1,nSy
     76            DO bi=1,nSx
     77               DO j=1,sNy
     78                  DO i=1,sNx
     79                     local(i,j,bi,bj) = xC(i,j,bi,bj)
     80                  ENDDO
     81               ENDDO
     82            ENDDO
     83         ENDDO
     84         CALL GATHER_2D_R8( xfer_array, local, myThid )
     85         _BEGIN_MASTER( myThid )
     86         buffsize = Nx*Ny
     87         CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
     88     &        0,10001005,toissmcomm,mpistatus)
     89         _END_MASTER( myThid )
     90         
     91C     Send latitude North of center of cell
     92         DO bj=1,nSy
     93            DO bi=1,nSx
     94               DO j=1,sNy
     95                  DO i=1,sNx
     96                     local(i,j,bi,bj) = yC(i,j,bi,bj)
     97                  ENDDO
     98               ENDDO
     99            ENDDO
     100         ENDDO
     101         CALL GATHER_2D_R8( xfer_array, local, myThid )
     102         _BEGIN_MASTER( myThid )
     103         buffsize = Nx*Ny
     104         CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
     105     &        0,10001006,toissmcomm,mpistatus)
     106         _END_MASTER( myThid )
     107
    80108      ENDIF
    81109       
Note: See TracChangeset for help on using the changeset viewer.