Changeset 22708
- Timestamp:
- 04/24/18 07:27:58 (7 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
TabularUnified issm/trunk-jpl/test/MITgcm/code/cpl_issm.F ¶
r22681 r22708 51 51 _R8 local(1:sNx,1:sNy,nSx,nSy) 52 52 53 C Initialization steps I1, I2, and I3: 53 54 IF( myTime .EQ. startTime ) THEN 54 55 55 C Send/receive scalar values56 C I1. ISSM sends CouplingTime, the interval at which we couple 56 57 IF( myProcId .EQ. 0 ) THEN 57 _BEGIN_MASTER( myThid ) 58 _BEGIN_MASTER( myThid ) 58 59 call MPI_Recv(CouplingTime,1,MPI_DOUBLE,0,10001000, 59 60 & toissmcomm,mpistatus,mpiRC) 60 61 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 65 C I2. MITgcm sends grid size (NX and NY) 66 IF( myProcId .EQ. 0 ) THEN 67 _BEGIN_MASTER( myThid ) 66 68 call MPI_Send(Nx,1,MPI_INT,0,10001003, 67 69 & toissmcomm,mpistatus) … … 69 71 & toissmcomm,mpistatus) 70 72 _END_MASTER( myThid ) 71 endif73 ENDIF 72 74 75 C I3. MITgcm sends grid coordinates of center of cells 76 C (longitude -180 <= XC < 180 and latitude YC) 73 77 C Send longitude East of center of cell 74 print*,'this is xC',xC75 78 DO bj=1,nSy 76 79 DO bi=1,nSx … … 82 85 ENDDO 83 86 ENDDO 84 print*,'this is local xC',local85 87 CALL BAR2( myThid ) 86 88 CALL GATHER_2D_R8( xfer_array, local, Nx, Ny, 87 89 & .FALSE., .FALSE., myThid ) 88 print*,'this is global XC',xfer_array89 90 IF( myProcId .EQ. 0 ) THEN 90 91 _BEGIN_MASTER( myThid ) … … 95 96 ENDIF 96 97 CALL BAR2( myThid ) 97 98 98 C Send latitude North of center of cell 99 print*,'this is yC',yC100 99 DO bj=1,nSy 101 100 DO bi=1,nSx … … 107 106 ENDDO 108 107 ENDDO 109 print*,'this is local yC',local110 108 CALL BAR2( myThid ) 111 109 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 ) 114 111 IF( myProcId .EQ. 0 ) THEN 115 112 _BEGIN_MASTER( myThid ) … … 120 117 ENDIF 121 118 CALL BAR2( myThid ) 122 print*,'Done Sending XC/YC arrays.'123 119 124 C Receive icebase 120 ENDIF 121 C End initialization steps I1, I2, and I3. 122 123 C Recurring steps C1 and C2: 124 IF( MOD(myTime,CouplingTime) .LT. deltaT/2. ) THEN 125 126 127 C C1. ISSM sends ice model time IceTimeTag 125 128 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_array129 _BEGIN_MASTER( myThid ) 130 call MPI_Recv(IceModelTime,1,MPI_DOUBLE,0,10001001, 131 & toissmcomm,mpistatus,mpiRC) 132 print*, 'Ocean received IceModelTime: ', IceModelTime 130 133 _END_MASTER( myThid ) 131 134 ENDIF 132 135 133 C Send melt rate 134 print*,'this is yC',yC 136 C 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 145 C End recurring steps C1 and C2. 146 147 C Recurring step C3 except during Initialization: 148 C C3. MITgcm sends 149 C (N-1)*CouplingTime <= OceanModelTime < N*CouplingTime 150 C time-mean melt rate to ISSM 151 IF( myTime .NE. startTime .AND. 152 & MOD(myTime,CouplingTime) .LT. deltaT/2. ) THEN 135 153 DO bj=1,nSy 136 154 DO bi=1,nSx … … 142 160 ENDDO 143 161 ENDDO 144 print*,'this is local shelficeHeatFlux',local145 162 CALL BAR2( myThid ) 146 163 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 ) 149 165 IF( myProcId .EQ. 0 ) THEN 150 166 _BEGIN_MASTER( myThid ) … … 156 172 CALL BAR2( myThid ) 157 173 print*,'Done Sending shelficeHeatFlux array.' 174 175 ENDIF 176 C End recurring step C3. 158 177 178 C Recurring step C4 except during Termination: 179 C 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 159 189 ENDIF 190 C End recurring step C4. 160 191 161 #endif /* ALLOW_CPL_ MPMICE*/192 #endif /* ALLOW_CPL_ISSM */ 162 193 163 194 RETURN
Note:
See TracChangeset
for help on using the changeset viewer.