/[MITgcm]/MITgcm/model/src/ini_depths.F
ViewVC logotype

Contents of /MITgcm/model/src/ini_depths.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.39 - (show annotations) (download)
Sun Nov 5 18:36:05 2006 UTC (17 years, 6 months ago) by edhill
Branch: MAIN
CVS Tags: checkpoint59, checkpoint58y_post, checkpoint58t_post, checkpoint58w_post, mitgcm_mapl_00, checkpoint58r_post, checkpoint59p, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint58v_post, checkpoint58x_post, checkpoint59j, checkpoint58u_post, checkpoint58s_post
Changes since 1.38: +24 -2 lines
add ability to read bathy, salt, and theta using MNC (off by def)
  -- and this can be readily extended to most of the other files
  in PARM05 of the main "data" namelist file

1 C $Header: /u/gcmpack/MITgcm/model/src/ini_depths.F,v 1.38 2006/10/17 18:52:34 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: INI_DEPTHS
9 C !INTERFACE:
10 SUBROUTINE INI_DEPTHS( myThid )
11 C !DESCRIPTION: \bv
12 C *==========================================================*
13 C | SUBROUTINE INI_DEPTHS
14 C | o define R_position of Lower and Surface Boundaries
15 C *==========================================================*
16 C |atmosphere orography:
17 C | define either in term of P_topo or converted from Z_topo
18 C |ocean bathymetry:
19 C | The depths of the bottom of the model is specified in
20 C | terms of an XY map with one depth for each column of
21 C | grid cells. Depths do not have to coincide with the
22 C | model levels. The model lopping algorithm makes it
23 C | possible to represent arbitrary depths.
24 C | The mode depths map also influences the models topology
25 C | By default the model domain wraps around in X and Y.
26 C | This default doubly periodic topology is "supressed"
27 C | if a depth map is defined which closes off all wrap
28 C | around flow.
29 C *==========================================================*
30 C \ev
31
32 C !USES:
33 IMPLICIT NONE
34 C === Global variables ===
35 #include "SIZE.h"
36 #include "EEPARAMS.h"
37 #include "PARAMS.h"
38 #include "GRID.h"
39 #include "SURFACE.h"
40 #ifdef ALLOW_SHELFICE
41 # include "SHELFICE.h"
42 #endif /* ALLOW_SHELFICE */
43 #ifdef ALLOW_MNC
44 #include "MNC_PARAMS.h"
45 #endif
46
47 C !INPUT/OUTPUT PARAMETERS:
48 C == Routine arguments ==
49 C myThid - Number of this instance of INI_DEPTHS
50 INTEGER myThid
51 CEndOfInterface
52
53 C !LOCAL VARIABLES:
54 C == Local variables ==
55 C iG, jG - Global coordinate index
56 C bi,bj - Tile indices
57 C I,J,K - Loop counters
58 C oldPrec - Temporary used in controlling binary input dataset precision
59 C msgBuf - Informational/error meesage buffer
60 INTEGER iG, jG
61 INTEGER bi, bj
62 INTEGER I, J
63 CHARACTER*(MAX_LEN_MBUF) msgBuf
64 CEOP
65
66 IF (usingPCoords .AND. bathyFile .NE. ' '
67 & .AND. topoFile .NE. ' ' ) THEN
68 WRITE(msgBuf,'(A,A)')
69 & 'S/R INI_DEPTHS: both bathyFile & topoFile are specified:',
70 & ' select the right one !'
71 CALL PRINT_ERROR( msgBuf , myThid)
72 STOP 'ABNORMAL END: S/R INI_DEPTHS'
73 ENDIF
74
75 C------
76 C 0) Initialize R_low and Ro_surf (define an empty domain)
77 C------
78 DO bj = myByLo(myThid), myByHi(myThid)
79 DO bi = myBxLo(myThid), myBxHi(myThid)
80 DO j=1-Oly,sNy+Oly
81 DO i=1-Olx,sNx+Olx
82 R_low(i,j,bi,bj) = 0. _d 0
83 Ro_surf(i,j,bi,bj) = 0. _d 0
84 topoZ(i,j,bi,bj) = 0. _d 0
85 #ifdef ALLOW_SHELFICE
86 R_shelfIce(i,j,bi,bj) = 0. _d 0
87 #endif /* ALLOW_SHELFICE */
88 ENDDO
89 ENDDO
90 ENDDO
91 ENDDO
92
93 C- Need to synchronize here before doing master-thread IO
94 _BARRIER
95
96 C------
97 C 1) Set R_low = the Lower (in r sense) boundary of the fluid column :
98 C------
99 IF (usingPCoords .OR. bathyFile .EQ. ' ') THEN
100 C- e.g., atmosphere : R_low = Top of atmosphere
101 C- ocean : R_low = Bottom
102 DO bj = myByLo(myThid), myByHi(myThid)
103 DO bi = myBxLo(myThid), myBxHi(myThid)
104 DO j=1,sNy
105 DO i=1,sNx
106 R_low(i,j,bi,bj) = rF(Nr+1)
107 ENDDO
108 ENDDO
109 ENDDO
110 ENDDO
111 ELSE
112
113 #ifdef ALLOW_MNC
114 IF (useMNC .AND. mnc_read_bathy) THEN
115 CALL MNC_CW_ADD_VNAME('bathy', 'Cen_xy_Hn__-__-', 3,4, myThid)
116 CALL MNC_FILE_CLOSE_ALL_MATCHING(bathyFile, myThid)
117 CALL MNC_CW_SET_UDIM(bathyFile, 1, myThid)
118 CALL MNC_CW_SET_CITER(bathyFile, 2, -1, -1, -1, myThid)
119 CALL MNC_CW_SET_UDIM(bathyFile, 1, myThid)
120 CALL MNC_CW_RL_R('D',bathyFile,0,0,'bathy',R_low, myThid)
121 CALL MNC_FILE_CLOSE_ALL_MATCHING(bathyFile, myThid)
122 CALL MNC_CW_DEL_VNAME('bathy', myThid)
123 ELSE
124 #endif /* ALLOW_MNC */
125 C Read the bathymetry using the mid-level I/O package read_write_rec
126 C The 0 is the "iteration" argument. The 1 is the record number.
127 CALL READ_REC_XY_RS( bathyFile, R_low, 1, 0, myThid )
128 C Read the bathymetry using the mid-level I/O package read_write_fld
129 C The 0 is the "iteration" argument. The ' ' is an empty suffix
130 c CALL READ_FLD_XY_RS( bathyFile, ' ', R_low, 0, myThid )
131 C Read the bathymetry using the low-level I/O package
132 c CALL MDSREADFIELD( bathyFile, readBinaryPrec,
133 c & 'RS', 1, R_low, 1, myThid )
134
135 #ifdef ALLOW_MNC
136 ENDIF
137 #endif /* ALLOW_MNC */
138
139
140 ENDIF
141 C- end setup R_low in the interior
142
143 C- fill in the overlap (+ BARRIER):
144 _EXCH_XY_R4(R_low, myThid )
145
146 IF ( debugLevel.GE.debLevB ) THEN
147 c PRINT *, ' Calling plot field', myThid
148 CALL PLOT_FIELD_XYRS( R_low, 'Bottom depths (ini_depths)',
149 & -1, myThid )
150 ENDIF
151 c CALL WRITE_FLD_XY_RS( 'R_low' ,' ', R_low, 0,myThid)
152
153 c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
154
155 C------
156 C 2) Set R_surf = Surface boundary: ocean surface / ground for the atmosphere
157 C------
158
159 IF ( usingPCoords .AND. bathyFile.NE.' ' ) THEN
160 C------ read directly Po_surf from bathyFile (only for backward compatibility)
161
162 CALL READ_REC_XY_RS( bathyFile, Ro_surf, 1, 0, myThid )
163
164 ELSEIF ( topoFile.EQ.' ' ) THEN
165 C------ set default value:
166
167 DO bj = myByLo(myThid), myByHi(myThid)
168 DO bi = myBxLo(myThid), myBxHi(myThid)
169 DO j=1,sNy
170 DO i=1,sNx
171 Ro_surf(i,j,bi,bj) = Ro_SeaLevel
172 ENDDO
173 ENDDO
174 ENDDO
175 ENDDO
176
177 ELSE
178 C------ read from file:
179
180 C- read surface topography (in m) from topoFile (case topoFile.NE.' '):
181
182 CALL READ_REC_XY_RS( topoFile, topoZ, 1, 0, myThid )
183 _BARRIER
184
185 IF (buoyancyRelation .EQ. 'ATMOSPHERIC') THEN
186 C----
187 C Convert Surface Geopotential to (reference) Surface Pressure
188 C according to Tref profile, using same discretisation as in calc_phi_hyd
189 C----
190 c CALL WRITE_FLD_XY_RS( 'topo_Z',' ',topoZ,0,myThid)
191
192 CALL INI_P_GROUND( 2, topoZ,
193 O Ro_surf,
194 I myThid )
195
196 c _BARRIER
197 C This I/O is now done in write_grid.F
198 c CALL WRITE_FLD_XY_RS( 'topo_P',' ',Ro_surf,0,myThid)
199
200 ELSE
201 C----
202 C Direct Transfer to Ro_surf (e.g., to specify upper ocean boundary
203 C below an ice-shelf - NOTE - actually not yet implemented )
204 DO bj = myByLo(myThid), myByHi(myThid)
205 DO bi = myBxLo(myThid), myBxHi(myThid)
206 DO j=1,sNy
207 DO i=1,sNx
208 Ro_surf(i,j,bi,bj) = topoZ(i,j,bi,bj)
209 ENDDO
210 ENDDO
211 ENDDO
212 ENDDO
213
214 ENDIF
215
216 C------ end case "read topoFile"
217 ENDIF
218
219 C----- fill in the overlap (+ BARRIER):
220 _EXCH_XY_R4(Ro_surf, myThid )
221
222 c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
223
224 C------
225 C 3) Close the Domain (special configuration).
226 C------
227 IF (usingPCoords) THEN
228 DO bj = myByLo(myThid), myByHi(myThid)
229 DO bi = myBxLo(myThid), myBxHi(myThid)
230 DO j=1-Oly,sNy+Oly
231 DO i=1-Olx,sNx+Olx
232 iG = myXGlobalLo-1+(bi-1)*sNx+I
233 jG = myYGlobalLo-1+(bj-1)*sNy+J
234 C Test for eastern edge
235 c IF ( iG .EQ. Nx ) Ro_surf(i,j,bi,bj) = 0.
236 C Test for northern edge
237 c IF ( jG .EQ. Ny ) Ro_surf(i,j,bi,bj) = 0.
238 IF (usingSphericalPolarGrid .AND. ABS(yC(I,J,bi,bj)).GE.90. )
239 & Ro_surf(I,J,bi,bj) = rF(Nr+1)
240 ENDDO
241 ENDDO
242 ENDDO
243 ENDDO
244 ELSE
245 DO bj = myByLo(myThid), myByHi(myThid)
246 DO bi = myBxLo(myThid), myBxHi(myThid)
247 DO j=1-Oly,sNy+Oly
248 DO i=1-Olx,sNx+Olx
249 iG = myXGlobalLo-1+(bi-1)*sNx+I
250 jG = myYGlobalLo-1+(bj-1)*sNy+J
251 C Test for eastern edge
252 c IF ( iG .EQ. Nx ) R_low(i,j,bi,bj) = 0.
253 C Test for northern edge
254 c IF ( jG .EQ. Ny ) R_low(i,j,bi,bj) = 0.
255 IF (usingSphericalPolarGrid .AND. ABS(yC(I,J,bi,bj)).GE.90. )
256 & R_low(I,J,bi,bj) = Ro_SeaLevel
257 ENDDO
258 ENDDO
259 ENDDO
260 ENDDO
261 ENDIF
262
263 IF ( debugLevel.GE.debLevB ) THEN
264 CALL PLOT_FIELD_XYRS( Ro_surf,
265 & 'Surface reference r-position (ini_depths)',
266 & -1, myThid )
267 ENDIF
268 c CALL WRITE_FLD_XY_RS('Ro_surf',' ',Ro_surf,0,myThid)
269
270 #ifdef ALLOW_SHELFICE
271 c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
272 IF ( useShelfIce ) THEN
273 C------
274 C 4) Set R_shelfIce = the Lower (in r sense) boundary of floating shelfice :
275 C------
276 IF (usingPCoords .OR. shelfIceFile .EQ. ' ') THEN
277 C- e.g., atmosphere : R_low = Top of atmosphere
278 C- ocean : R_low = Bottom
279 DO bj = myByLo(myThid), myByHi(myThid)
280 DO bi = myBxLo(myThid), myBxHi(myThid)
281 DO j=1,sNy
282 DO i=1,sNx
283 R_shelfIce(i,j,bi,bj) = 0. _d 0
284 ENDDO
285 ENDDO
286 ENDDO
287 ENDDO
288 ELSE
289 C Read the shelfIce draught using the mid-level I/O pacakage read_write_rec
290 C The 0 is the "iteration" argument. The 1 is the record number.
291 CALL READ_REC_XY_RS( shelfIceFile, R_shelfIce, 1, 0, myThid )
292 C Read the shelfIce draught using the mid-level I/O pacakage read_write_fld
293 C The 0 is the "iteration" argument. The ' ' is an empty suffix
294 C CALL READ_FLD_XY_RS( shelfIceFile, ' ', R_shelfIce, 0, myThid )
295 c CALL READ_FLD_XY_RS( bathyFile, ' ', R_low, 0, myThid )
296 C Read the selfIce draught using the low-level I/O package
297 c CALL MDSREADFIELD( shelfIceFile, readBinaryPrec,
298 c & 'RS', 1, R_selfIce, 1, myThid )
299
300 ENDIF
301 C- end setup R_shelfIce in the interior
302
303 C- fill in the overlap (+ BARRIER):
304 _EXCH_XY_R4(R_shelfIce, myThid )
305
306 c CALL PLOT_FIELD_XYRS(R_selfIce,'Shelf ice draught (ini_depths)',
307 c & 1,myThid)
308 CML CALL WRITE_FLD_XY_RS( 'R_shelfIce' ,' ', R_shelfIce, 0,myThid)
309 ENDIF
310 #endif /* ALLOW_SHELFICE */
311
312 C-- Everyone else must wait for the depth to be loaded
313 C- note: might not be necessary since all single-thread IO above
314 C are followed by an EXCH (with BARRIER)
315 _BARRIER
316
317 RETURN
318 END

  ViewVC Help
Powered by ViewVC 1.1.22