/[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.44 - (show annotations) (download)
Fri Apr 24 01:55:03 2009 UTC (15 years, 1 month ago) by jmc
Branch: MAIN
Changes since 1.43: +7 -2 lines
rename OBCS_CHECK_TOPOGRAPHY -> OBCS_CHECK_DEPTHS
+ call EXCH2_CHECK_DEPTHS

1 C $Header: /u/gcmpack/MITgcm/model/src/ini_depths.F,v 1.43 2008/10/02 00:33:20 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_MNC
41 # include "MNC_PARAMS.h"
42 #endif
43
44 C !INPUT/OUTPUT PARAMETERS:
45 C == Routine arguments ==
46 C myThid - Number of this instance of INI_DEPTHS
47 INTEGER myThid
48 CEndOfInterface
49
50 C !LOCAL VARIABLES:
51 C == Local variables ==
52 C iG, jG - Global coordinate index
53 C bi,bj - Tile indices
54 C I,J,K - Loop counters
55 C oldPrec - Temporary used in controlling binary input dataset precision
56 C msgBuf - Informational/error meesage buffer
57 INTEGER iG, jG
58 INTEGER bi, bj
59 INTEGER I, J
60 CHARACTER*(MAX_LEN_MBUF) msgBuf
61 CEOP
62
63 IF (usingPCoords .AND. bathyFile .NE. ' '
64 & .AND. topoFile .NE. ' ' ) THEN
65 WRITE(msgBuf,'(A,A)')
66 & 'S/R INI_DEPTHS: both bathyFile & topoFile are specified:',
67 & ' select the right one !'
68 CALL PRINT_ERROR( msgBuf , myThid)
69 STOP 'ABNORMAL END: S/R INI_DEPTHS'
70 ENDIF
71
72 C------
73 C 0) Initialize R_low and Ro_surf (define an empty domain)
74 C------
75 DO bj = myByLo(myThid), myByHi(myThid)
76 DO bi = myBxLo(myThid), myBxHi(myThid)
77 DO j=1-Oly,sNy+Oly
78 DO i=1-Olx,sNx+Olx
79 R_low(i,j,bi,bj) = 0. _d 0
80 Ro_surf(i,j,bi,bj) = 0. _d 0
81 topoZ(i,j,bi,bj) = 0. _d 0
82 ENDDO
83 ENDDO
84 ENDDO
85 ENDDO
86
87 C- Need to synchronize here before doing master-thread IO
88 _BARRIER
89
90 C------
91 C 1) Set R_low = the Lower (in r sense) boundary of the fluid column :
92 C------
93 IF (usingPCoords .OR. bathyFile .EQ. ' ') THEN
94 C- e.g., atmosphere : R_low = Top of atmosphere
95 C- ocean : R_low = Bottom
96 DO bj = myByLo(myThid), myByHi(myThid)
97 DO bi = myBxLo(myThid), myBxHi(myThid)
98 DO j=1,sNy
99 DO i=1,sNx
100 R_low(i,j,bi,bj) = rF(Nr+1)
101 ENDDO
102 ENDDO
103 ENDDO
104 ENDDO
105 ELSE
106
107 #ifdef ALLOW_MNC
108 IF (useMNC .AND. mnc_read_bathy) THEN
109 CALL MNC_CW_ADD_VNAME('bathy', 'Cen_xy_Hn__-__-', 3,4, myThid)
110 CALL MNC_FILE_CLOSE_ALL_MATCHING(bathyFile, myThid)
111 CALL MNC_CW_SET_UDIM(bathyFile, 1, myThid)
112 CALL MNC_CW_SET_CITER(bathyFile, 2, -1, -1, -1, myThid)
113 CALL MNC_CW_SET_UDIM(bathyFile, 1, myThid)
114 CALL MNC_CW_RL_R('D',bathyFile,0,0,'bathy',R_low, myThid)
115 CALL MNC_FILE_CLOSE_ALL_MATCHING(bathyFile, myThid)
116 CALL MNC_CW_DEL_VNAME('bathy', myThid)
117 ELSE
118 #endif /* ALLOW_MNC */
119 C Read the bathymetry using the mid-level I/O package read_write_rec
120 C The 0 is the "iteration" argument. The 1 is the record number.
121 CALL READ_REC_XY_RS( bathyFile, R_low, 1, 0, myThid )
122 C Read the bathymetry using the mid-level I/O package read_write_fld
123 C The 0 is the "iteration" argument. The ' ' is an empty suffix
124 c CALL READ_FLD_XY_RS( bathyFile, ' ', R_low, 0, myThid )
125 C Read the bathymetry using the low-level I/O package
126 c CALL MDSREADFIELD( bathyFile, readBinaryPrec,
127 c & 'RS', 1, R_low, 1, myThid )
128
129 #ifdef ALLOW_MNC
130 ENDIF
131 #endif /* ALLOW_MNC */
132
133
134 ENDIF
135 C- end setup R_low in the interior
136
137 C- fill in the overlap (+ BARRIER):
138 _EXCH_XY_R4(R_low, myThid )
139
140 IF ( debugLevel.GE.debLevB ) THEN
141 c PRINT *, ' Calling plot field', myThid
142 CALL PLOT_FIELD_XYRS( R_low, 'Bottom depths (ini_depths)',
143 & -1, myThid )
144 ENDIF
145 c CALL WRITE_FLD_XY_RS( 'R_low' ,' ', R_low, 0,myThid)
146
147 c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
148
149 C------
150 C 2) Set R_surf = Surface boundary: ocean surface / ground for the atmosphere
151 C------
152
153 IF ( usingPCoords .AND. bathyFile.NE.' ' ) THEN
154 C------ read directly Po_surf from bathyFile (only for backward compatibility)
155
156 CALL READ_REC_XY_RS( bathyFile, Ro_surf, 1, 0, myThid )
157
158 ELSEIF ( topoFile.EQ.' ' ) THEN
159 C------ set default value:
160
161 DO bj = myByLo(myThid), myByHi(myThid)
162 DO bi = myBxLo(myThid), myBxHi(myThid)
163 DO j=1,sNy
164 DO i=1,sNx
165 Ro_surf(i,j,bi,bj) = rF(1)
166 ENDDO
167 ENDDO
168 ENDDO
169 ENDDO
170
171 ELSE
172 C------ read from file:
173
174 C- read surface topography (in m) from topoFile (case topoFile.NE.' '):
175
176 CALL READ_REC_XY_RS( topoFile, topoZ, 1, 0, myThid )
177 _BARRIER
178
179 IF (buoyancyRelation .EQ. 'ATMOSPHERIC') THEN
180 C----
181 C Convert Surface Geopotential to (reference) Surface Pressure
182 C according to Tref profile, using same discretisation as in calc_phi_hyd
183 C----
184 c CALL WRITE_FLD_XY_RS( 'topo_Z',' ',topoZ,0,myThid)
185
186 CALL INI_P_GROUND( 2, topoZ,
187 O Ro_surf,
188 I myThid )
189
190 c _BARRIER
191 C This I/O is now done in write_grid.F
192 c CALL WRITE_FLD_XY_RS( 'topo_P',' ',Ro_surf,0,myThid)
193
194 ELSEIF ( buoyancyRelation.EQ.'OCEANICP' ) THEN
195
196 WRITE(msgBuf,'(A,A)') 'S/R INI_DEPTHS: ',
197 & 'from topoFile (in m) to ref.bottom pressure: Not yet coded'
198 CALL PRINT_ERROR( msgBuf , myThid)
199 STOP 'ABNORMAL END: S/R INI_DEPTHS'
200
201 ELSE
202 C----
203 C Direct Transfer to Ro_surf (e.g., to specify upper ocean boundary
204 C below an ice-shelf - NOTE - actually not yet implemented )
205 DO bj = myByLo(myThid), myByHi(myThid)
206 DO bi = myBxLo(myThid), myBxHi(myThid)
207 DO j=1,sNy
208 DO i=1,sNx
209 Ro_surf(i,j,bi,bj) = topoZ(i,j,bi,bj)
210 ENDDO
211 ENDDO
212 ENDDO
213 ENDDO
214
215 ENDIF
216
217 C------ end case "read topoFile"
218 ENDIF
219
220 C----- fill in the overlap (+ BARRIER):
221 _EXCH_XY_R4(Ro_surf, myThid )
222
223 c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
224
225 C------
226 C 3) Close the Domain (special configuration).
227 C------
228 IF (usingPCoords) THEN
229 DO bj = myByLo(myThid), myByHi(myThid)
230 DO bi = myBxLo(myThid), myBxHi(myThid)
231 DO j=1-Oly,sNy+Oly
232 DO i=1-Olx,sNx+Olx
233 iG = myXGlobalLo-1+(bi-1)*sNx+I
234 jG = myYGlobalLo-1+(bj-1)*sNy+J
235 C Test for eastern edge
236 c IF ( iG .EQ. Nx ) Ro_surf(i,j,bi,bj) = 0.
237 C Test for northern edge
238 c IF ( jG .EQ. Ny ) Ro_surf(i,j,bi,bj) = 0.
239 IF (usingSphericalPolarGrid .AND. ABS(yC(I,J,bi,bj)).GE.90. )
240 & Ro_surf(I,J,bi,bj) = rF(Nr+1)
241 ENDDO
242 ENDDO
243 ENDDO
244 ENDDO
245 ELSE
246 DO bj = myByLo(myThid), myByHi(myThid)
247 DO bi = myBxLo(myThid), myBxHi(myThid)
248 DO j=1-Oly,sNy+Oly
249 DO i=1-Olx,sNx+Olx
250 iG = myXGlobalLo-1+(bi-1)*sNx+I
251 jG = myYGlobalLo-1+(bj-1)*sNy+J
252 C Test for eastern edge
253 c IF ( iG .EQ. Nx ) R_low(i,j,bi,bj) = 0.
254 C Test for northern edge
255 c IF ( jG .EQ. Ny ) R_low(i,j,bi,bj) = 0.
256 IF (usingSphericalPolarGrid .AND. ABS(yC(I,J,bi,bj)).GE.90. )
257 & R_low(I,J,bi,bj) = rF(1)
258 ENDDO
259 ENDDO
260 ENDDO
261 ENDDO
262 ENDIF
263
264 IF ( debugLevel.GE.debLevB ) THEN
265 CALL PLOT_FIELD_XYRS( Ro_surf,
266 & 'Surface reference r-position (ini_depths)',
267 & -1, myThid )
268 ENDIF
269 c CALL WRITE_FLD_XY_RS('Ro_surf',' ',Ro_surf,0,myThid)
270
271 C-- Everyone else must wait for the depth to be loaded
272 C- note: might not be necessary since all single-thread IO above
273 C are followed by an EXCH (with BARRIER)
274 _BARRIER
275
276 #ifdef ALLOW_OBCS
277 IF ( useOBCS ) THEN
278 C check for inconsistent topography along boundaries and fix it
279 CALL OBCS_CHECK_DEPTHS( myThid )
280 C update the overlaps
281 _EXCH_XY_R4(R_low, myThid )
282 ENDIF
283 #endif /* ALLOW_OBCS */
284
285 #ifdef ALLOW_EXCH2
286 C Check domain boundary (e.g., in case of missing tiles)
287 CALL EXCH2_CHECK_DEPTHS( R_low, Ro_surf, myThid )
288 #endif
289
290 RETURN
291 END

  ViewVC Help
Powered by ViewVC 1.1.22