source: issm/trunk/test/MITgcm/code/cpl_issm.F@ 23189

Last change on this file since 23189 was 23189, checked in by Mathieu Morlighem, 7 years ago

merged trunk-jpl and trunk for revision 23187

File size: 7.0 KB
Line 
1#include "PACKAGES_CONFIG.h"
2#include "CPP_OPTIONS.h"
3
4CBOP
5C !ROUTINE: CPL_ISSM
6C !INTERFACE:
7 SUBROUTINE CPL_ISSM( myTime, myIter, myThid )
8
9C !DESCRIPTION: \bv
10C *==================================================================
11C | SUBROUTINE cpl_issm
12C | o Couple MITgcm ocean model with ISSM ice sheet model
13C *==================================================================
14C \ev
15
16C !USES:
17 IMPLICIT NONE
18C == Global variables ==
19#include "SIZE.h"
20#include "EEPARAMS.h"
21#include "PARAMS.h"
22#include "DYNVARS.h"
23#include "GRID.h"
24#include "FFIELDS.h"
25#include "SHELFICE_OPTIONS.h"
26#include "SHELFICE.h"
27#ifdef ALLOW_EXF
28# include "EXF_OPTIONS.h"
29# include "EXF_FIELDS.h"
30#endif
31
32 LOGICAL DIFFERENT_MULTIPLE
33 EXTERNAL DIFFERENT_MULTIPLE
34
35C !LOCAL VARIABLES:
36C mytime - time counter for this thread (seconds)
37C myiter - iteration counter for this thread
38C mythid - thread number for this instance of the routine.
39 _RL mytime
40 INTEGER myiter, mythid
41CEOP
42
43#ifdef ALLOW_CPL_ISSM
44#include "EESUPPORT.h"
45 COMMON /CPL_MPI_ID/ mpiMyWid, toissmcomm
46 INTEGER mpiMyWid, toissmcomm, mpiRC
47 INTEGER mpistatus(MPI_STATUS_SIZE)
48 INTEGER i, j, bi, bj, buffsize
49 COMMON /CPL_ISSM_TIME/ CouplingTime
50 _R8 CouplingTime, IceModelTime
51 _R8 xfer_array(Nx,Ny)
52 _R8 local(1:sNx,1:sNy,nSx,nSy)
53 CHARACTER*(MAX_LEN_MBUF) suff
54
55C Initialization steps I1, I2, and I3:
56 IF( myTime .EQ. startTime ) THEN
57
58C I1. ISSM sends CouplingTime, the interval at which we couple
59 IF( myProcId .EQ. 0 ) THEN
60 _BEGIN_MASTER( myThid )
61 call MPI_Recv(CouplingTime,1,MPI_DOUBLE,0,10001000,
62 & toissmcomm,mpistatus,mpiRC)
63 _END_MASTER( myThid )
64 ENDIF
65 _BEGIN_MASTER( myThid )
66 CALL MPI_BCAST(CouplingTime,1,MPI_DOUBLE,0,
67 & MPI_COMM_MODEL,mpiRC)
68 _END_MASTER( myThid )
69C print*, 'Ocean received CouplingTime: ', CouplingTime
70
71C I2. MITgcm sends grid size (NX and NY)
72 IF( myProcId .EQ. 0 ) THEN
73 _BEGIN_MASTER( myThid )
74 call MPI_Send(Nx,1,MPI_INT,0,10001003,
75 & toissmcomm,mpistatus)
76 call MPI_Send(Ny,1,MPI_INT,0,10001004,
77 & toissmcomm,mpistatus)
78 _END_MASTER( myThid )
79 ENDIF
80
81C I3. MITgcm sends grid coordinates of center of cells
82C (longitude -180 <= XC < 180 and latitude YC)
83C Send longitude East of center of cell
84 DO bj=1,nSy
85 DO bi=1,nSx
86 DO j=1,sNy
87 DO i=1,sNx
88 local(i,j,bi,bj) = xC(i,j,bi,bj)
89 ENDDO
90 ENDDO
91 ENDDO
92 ENDDO
93 CALL BAR2( myThid )
94 CALL GATHER_2D_R8( xfer_array, local, Nx, Ny,
95 & .FALSE., .FALSE., myThid )
96 IF( myProcId .EQ. 0 ) THEN
97 _BEGIN_MASTER( myThid )
98 buffsize = Nx*Ny
99 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
100 & 0,10001005,toissmcomm,mpistatus)
101 _END_MASTER( myThid )
102 ENDIF
103 CALL BAR2( myThid )
104C Send latitude North of center of cell
105 DO bj=1,nSy
106 DO bi=1,nSx
107 DO j=1,sNy
108 DO i=1,sNx
109 local(i,j,bi,bj) = yC(i,j,bi,bj)
110 ENDDO
111 ENDDO
112 ENDDO
113 ENDDO
114 CALL BAR2( myThid )
115 CALL GATHER_2D_R8( xfer_array, local, Nx, Ny,
116 & .FALSE., .FALSE., myThid )
117 IF( myProcId .EQ. 0 ) THEN
118 _BEGIN_MASTER( myThid )
119 buffsize = Nx*Ny
120 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
121 & 0,10001006,toissmcomm,mpistatus)
122 _END_MASTER( myThid )
123 ENDIF
124 CALL BAR2( myThid )
125
126 ENDIF
127C End initialization steps I1, I2, and I3.
128
129C Recurring steps C1 and C2:
130 IF( MOD(myTime,CouplingTime) .LT. deltaT/2. ) THEN
131
132C C1. ISSM sends ice model time IceTimeTag
133 IF( myProcId .EQ. 0 ) THEN
134 _BEGIN_MASTER( myThid )
135 call MPI_Recv(IceModelTime,1,MPI_DOUBLE,0,10001001,
136 & toissmcomm,mpistatus,mpiRC)
137C print*, 'Ocean received IceModelTime: ', IceModelTime
138 _END_MASTER( myThid )
139 ENDIF
140
141C C2. MITgcm sends ocean model time OceanTimeTag
142 IF( myProcId .EQ. 0 ) THEN
143 _BEGIN_MASTER( myThid )
144 call MPI_Send(myTime,1,MPI_DOUBLE,0,10001002,
145 & toissmcomm,mpistatus)
146 _END_MASTER( myThid )
147 ENDIF
148
149 ENDIF
150C End recurring steps C1 and C2.
151
152C Recurring step C3 except during Initialization:
153C C3. MITgcm sends
154C (N-1)*CouplingTime <= OceanModelTime < N*CouplingTime
155C time-mean melt rate to ISSM
156 IF( myTime .NE. startTime .AND.
157 & MOD(myTime,CouplingTime) .LT. deltaT/2. ) THEN
158 DO bj=1,nSy
159 DO bi=1,nSx
160 DO j=1,sNy
161 DO i=1,sNx
162 local(i,j,bi,bj)=shelficeFreshWaterFlux(i,j,bi,bj)
163 ENDDO
164 ENDDO
165 ENDDO
166 ENDDO
167 CALL BAR2( myThid )
168 CALL GATHER_2D_R8( xfer_array, local, Nx, Ny,
169 & .FALSE., .FALSE., myThid )
170 IF( myProcId .EQ. 0 ) THEN
171 _BEGIN_MASTER( myThid )
172 buffsize = Nx*Ny
173 CALL MPI_SEND(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
174 & 0,10001007,toissmcomm,mpistatus)
175 _END_MASTER( myThid )
176 ENDIF
177 CALL BAR2( myThid )
178C print*,'Done Sending shelficeFreshWaterFlux array.'
179
180 ENDIF
181C End recurring step C3.
182
183C Recurring step C4 except during Termination:
184C C4. ISSM sends IceModelTime=(N-1)*CouplingTime base to MITgcm
185 IF( myTime .NE. endtime .AND.
186 & MOD(myTime,CouplingTime) .LT. deltaT/2. ) THEN
187 WRITE(suff,'(I10.10)') myIter
188 CALL WRITE_FLD_XY_RS( 'R_shelfIce1_',suff,R_shelfIce,-1,myThid)
189 IF( myProcId .EQ. 0 ) THEN
190 _BEGIN_MASTER( myThid )
191 call MPI_Recv(xfer_array,buffsize,MPI_DOUBLE_PRECISION,
192 & 0,10001008,toissmcomm,mpistatus,mpiRC)
193 _END_MASTER( myThid )
194 ENDIF
195 CALL BAR2( myThid )
196 CALL SCATTER_2D_R8( xfer_array, local, Nx, Ny,
197 & .FALSE., .FALSE., myThid )
198 DO bj = myByLo(myThid), myByHi(myThid)
199 DO bi = myBxLo(myThid), myBxHi(myThid)
200 DO j=1,sNy
201 DO i=1,sNx
202 IF( local(i,j,bi,bj).LT.9998 ) THEN
203 R_shelfIce(i,j,bi,bj) = local(i,j,bi,bj)
204 ELSE
205 R_shelfIce(i,j,bi,bj) = 0. _d 0
206 ENDIF
207 ENDDO
208 ENDDO
209 ENDDO
210 ENDDO
211C- fill in the overlap (+ BARRIER):
212 _EXCH_XY_RS( R_shelfIce, myThid )
213 CALL WRITE_FLD_XY_RS( 'R_shelfIce2_',suff,R_shelfIce,-1,myThid)
214 ENDIF
215C End recurring step C4.
216
217#endif /* ALLOW_CPL_ISSM */
218
219 RETURN
220 END
Note: See TracBrowser for help on using the repository browser.