Changeset 23070


Ignore:
Timestamp:
08/07/18 12:18:38 (7 years ago)
Author:
dmenemen
Message:

MITgcm uses "base" sent by ISSM
code includes some diagnostics write statements, which will eventually
be commented out, and a manual check and removal of 1-m ice shelf
thickness, which will eventually be replaced to check for NaN value.
Verification output changes because there are truncation differences
at the <1e-5 level.

File:
1 edited

Legend:

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

    r22797 r23070  
    6666     &        MPI_COMM_MODEL,mpiRC)
    6767         _END_MASTER( myThid )
    68          print*, 'Ocean received CouplingTime: ', CouplingTime
     68C        print*, 'Ocean received CouplingTime: ', CouplingTime
    6969
    7070C   I2. MITgcm sends grid size (NX and NY)
     
    134134            call MPI_Recv(IceModelTime,1,MPI_DOUBLE,0,10001001,
    135135     &           toissmcomm,mpistatus,mpiRC)
    136             print*, 'Ocean received IceModelTime: ', IceModelTime
     136C           print*, 'Ocean received IceModelTime: ', IceModelTime
    137137            _END_MASTER( myThid )
    138138         ENDIF
     
    175175         ENDIF
    176176         CALL BAR2( myThid )
    177          print*,'Done Sending shelficeFreshWaterFlux array.'
     177C        print*,'Done Sending shelficeFreshWaterFlux array.'
    178178         
    179179      ENDIF
     
    182182C Recurring step C4 except during Termination:
    183183C  C4. ISSM sends IceModelTime=(N-1)*CouplingTime base to MITgcm
    184       IF( myTime .NE. endtime .AND.
    185      &     MOD(myTime,CouplingTime) .LT. deltaT/2. ) THEN         
     184CDM   IF( myTime .NE. endtime .AND.
     185      IF( myTime .EQ. startTime .AND.
     186     &     MOD(myTime,CouplingTime) .LT. deltaT/2. ) THEN
     187         CALL WRITE_FLD_XY_RS( 'R_shelfIce_1',' ',R_shelfIce,-1,myThid)
    186188         IF( myProcId .EQ. 0 ) THEN
    187189            _BEGIN_MASTER( myThid )         
     
    190192            _END_MASTER( myThid )
    191193         ENDIF
     194         CALL BAR2( myThid )
     195         CALL SCATTER_2D_R8( xfer_array, local, Nx, Ny,
     196     &        .FALSE., .FALSE., myThid )
     197         DO bj = myByLo(myThid), myByHi(myThid)
     198            DO bi = myBxLo(myThid), myBxHi(myThid)
     199               DO j=1,sNy
     200                  DO i=1,sNx
     201                     IF( ABS(local(i,j,bi,bj)+
     202     &                    0.89190188 ).GT.1E-5 ) THEN
     203                        R_shelfIce(i,j,bi,bj) = local(i,j,bi,bj)
     204                     ENDIF
     205                  ENDDO
     206               ENDDO
     207            ENDDO
     208         ENDDO
     209C- fill in the overlap (+ BARRIER):
     210         _EXCH_XY_RS( R_shelfIce, myThid )
     211         CALL WRITE_FLD_XY_RS( 'R_shelfIce_2',' ',R_shelfIce,-1,myThid)
    192212      ENDIF
    193213C End recurring step C4.
Note: See TracChangeset for help on using the changeset viewer.