Changeset 22708


Ignore:
Timestamp:
04/24/18 07:27:58 (7 years ago)
Author:
dmenemen
Message:

implementing coupling rules

File:
1 edited

Legend:

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

    r22681 r22708  
    5151      _R8 local(1:sNx,1:sNy,nSx,nSy)
    5252
     53C Initialization steps I1, I2, and I3:
    5354      IF( myTime .EQ. startTime ) THEN
    5455
    55 C     Send/receive scalar values
     56C   I1. ISSM sends CouplingTime, the interval at which we couple
    5657         IF( myProcId .EQ. 0 ) THEN
    57             _BEGIN_MASTER( myThid )         
     58            _BEGIN_MASTER( myThid )
    5859            call MPI_Recv(CouplingTime,1,MPI_DOUBLE,0,10001000,
    5960     &           toissmcomm,mpistatus,mpiRC)
    6061            print*, 'Ocean received CouplingTime: ', CouplingTime
    61             call MPI_Recv(IceModelTime,1,MPI_DOUBLE,0,10001001,
    62      &           toissmcomm,mpistatus,mpiRC)
    63             print*, 'Ocean received IceModelTime: ', IceModelTime
    64             call MPI_Send(myTime,1,MPI_DOUBLE,0,10001002,
    65      &           toissmcomm,mpistatus)
     62            _END_MASTER( myThid )
     63         ENDIF
     64
     65C   I2. MITgcm sends grid size (NX and NY)
     66         IF( myProcId .EQ. 0 ) THEN
     67            _BEGIN_MASTER( myThid )
    6668            call MPI_Send(Nx,1,MPI_INT,0,10001003,
    6769     &           toissmcomm,mpistatus)
     
    6971     &           toissmcomm,mpistatus)
    7072            _END_MASTER( myThid )
    71          endif
     73         ENDIF
    7274
     75C   I3. MITgcm sends grid coordinates of center of cells
     76C       (longitude -180 <= XC < 180 and latitude YC)
    7377C     Send longitude East of center of cell
    74          print*,'this is xC',xC
    7578         DO bj=1,nSy
    7679            DO bi=1,nSx
     
    8285            ENDDO
    8386         ENDDO
    84          print*,'this is local xC',local
    8587         CALL BAR2( myThid )
    8688         CALL GATHER_2D_R8( xfer_array, local, Nx, Ny,
    8789     &        .FALSE., .FALSE., myThid )
    88          print*,'this is global XC',xfer_array
    8990         IF( myProcId .EQ. 0 ) THEN
    9091            _BEGIN_MASTER( myThid )
     
    9596         ENDIF
    9697         CALL BAR2( myThid )
    97 
    9898C     Send latitude North of center of cell
    99          print*,'this is yC',yC
    10099         DO bj=1,nSy
    101100            DO bi=1,nSx
     
    107106            ENDDO
    108107         ENDDO
    109          print*,'this is local yC',local
    110108         CALL BAR2( myThid )
    111109         CALL GATHER_2D_R8( xfer_array, local, Nx, Ny,
    112      &                      .FALSE., .FALSE., myThid )
    113          print*,'this is global YC',xfer_array
     110     &        .FALSE., .FALSE., myThid )
    114111         IF( myProcId .EQ. 0 ) THEN
    115112            _BEGIN_MASTER( myThid )
     
    120117         ENDIF
    121118         CALL BAR2( myThid )
    122          print*,'Done Sending XC/YC arrays.'
    123119
    124 C     Receive icebase
     120      ENDIF
     121C End initialization steps I1, I2, and I3.
     122
     123C Recurring steps C1 and C2:
     124      IF( MOD(myTime,CouplingTime) .LT. deltaT/2. ) THEN
     125
     126
     127C   C1. ISSM sends ice model time IceTimeTag
    125128         IF( myProcId .EQ. 0 ) THEN
    126             _BEGIN_MASTER( myThid )         
    127             call MPI_Recv(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
    128      &           0,10001008,toissmcomm,mpistatus,mpiRC)
    129             print*, 'Ocean received icebase',xfer_array
     129            _BEGIN_MASTER( myThid )
     130            call MPI_Recv(IceModelTime,1,MPI_DOUBLE,0,10001001,
     131     &           toissmcomm,mpistatus,mpiRC)
     132            print*, 'Ocean received IceModelTime: ', IceModelTime
    130133            _END_MASTER( myThid )
    131134         ENDIF
    132135
    133 C     Send melt rate
    134          print*,'this is yC',yC
     136C   C2. MITgcm sends ocean model time OceanTimeTag
     137         IF( myProcId .EQ. 0 ) THEN
     138            _BEGIN_MASTER( myThid )
     139            call MPI_Send(myTime,1,MPI_DOUBLE,0,10001002,
     140     &           toissmcomm,mpistatus)
     141            _END_MASTER( myThid )
     142         ENDIF
     143
     144      ENDIF
     145C End recurring steps C1 and C2.
     146
     147C Recurring step C3 except during Initialization:
     148C  C3. MITgcm sends
     149C      (N-1)*CouplingTime <= OceanModelTime < N*CouplingTime
     150C      time-mean melt rate to ISSM
     151      IF( myTime .NE. startTime .AND.
     152     &     MOD(myTime,CouplingTime) .LT. deltaT/2. ) THEN
    135153         DO bj=1,nSy
    136154            DO bi=1,nSx
     
    142160            ENDDO
    143161         ENDDO
    144          print*,'this is local shelficeHeatFlux',local
    145162         CALL BAR2( myThid )
    146163         CALL GATHER_2D_R8( xfer_array, local, Nx, Ny,
    147      &                      .FALSE., .FALSE., myThid )
    148          print*,'this is global shelficeHeatFlux',xfer_array
     164     &        .FALSE., .FALSE., myThid )
    149165         IF( myProcId .EQ. 0 ) THEN
    150166            _BEGIN_MASTER( myThid )
     
    156172         CALL BAR2( myThid )
    157173         print*,'Done Sending shelficeHeatFlux array.'
     174         
     175      ENDIF
     176C End recurring step C3.
    158177
     178C Recurring step C4 except during Termination:
     179C  C4. ISSM sends IceModelTime=(N-1)*CouplingTime base to MITgcm
     180      IF( myTime .NE. endtime .AND.
     181     &     MOD(myTime,CouplingTime) .LT. deltaT/2. ) THEN         
     182         IF( myProcId .EQ. 0 ) THEN
     183            _BEGIN_MASTER( myThid )         
     184            call MPI_Recv(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
     185     &           0,10001008,toissmcomm,mpistatus,mpiRC)
     186            print*, 'Ocean received icebase',xfer_array
     187            _END_MASTER( myThid )
     188         ENDIF
    159189      ENDIF
     190C End recurring step C4.
    160191
    161 #endif /* ALLOW_CPL_MPMICE */
     192#endif /* ALLOW_CPL_ISSM */
    162193
    163194      RETURN
Note: See TracChangeset for help on using the changeset viewer.