/[MITgcm]/MITgcm_contrib/cam_devel/sigma_testing/code-sigma/ini_depths.F
ViewVC logotype

Contents of /MITgcm_contrib/cam_devel/sigma_testing/code-sigma/ini_depths.F

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


Revision 1.1 - (show annotations) (download)
Wed Jan 6 04:31:15 2010 UTC (15 years, 7 months ago) by cnh
Branch: MAIN
CVS Tags: HEAD
start some cam work

1 C $Header: /u/gcmpack/MITgcm/model/src/ini_depths.F,v 1.47 2009/08/28 19:47:10 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 C this is done within IO routines => no longer needed
89 c _BARRIER
90
91 C------
92 C 1) Set R_low = the Lower (in r sense) boundary of the fluid column :
93 C------
94 IF (usingPCoords .OR. bathyFile .EQ. ' ') THEN
95 C- e.g., atmosphere : R_low = Top of atmosphere
96 C- ocean : R_low = Bottom
97 DO bj = myByLo(myThid), myByHi(myThid)
98 DO bi = myBxLo(myThid), myBxHi(myThid)
99 DO j=1,sNy
100 DO i=1,sNx
101 R_low(i,j,bi,bj) = rF(Nr+1)
102 ENDDO
103 ENDDO
104 ENDDO
105 ENDDO
106 ELSE
107
108 #ifdef ALLOW_MNC
109 IF (useMNC .AND. mnc_read_bathy) THEN
110 CALL MNC_CW_ADD_VNAME('bathy', 'Cen_xy_Hn__-__-', 3,4, myThid)
111 CALL MNC_FILE_CLOSE_ALL_MATCHING(bathyFile, myThid)
112 CALL MNC_CW_SET_UDIM(bathyFile, 1, myThid)
113 CALL MNC_CW_SET_CITER(bathyFile, 2, -1, -1, -1, myThid)
114 CALL MNC_CW_SET_UDIM(bathyFile, 1, myThid)
115 CALL MNC_CW_RS_R('D',bathyFile,0,0,'bathy',R_low, myThid)
116 CALL MNC_FILE_CLOSE_ALL_MATCHING(bathyFile, myThid)
117 CALL MNC_CW_DEL_VNAME('bathy', myThid)
118 ELSE
119 #endif /* ALLOW_MNC */
120 C Read the bathymetry using the mid-level I/O package read_write_rec
121 C The 0 is the "iteration" argument. The 1 is the record number.
122 CALL READ_REC_XY_RS( bathyFile, R_low, 1, 0, myThid )
123 C Read the bathymetry using the mid-level I/O package read_write_fld
124 C The 0 is the "iteration" argument. The ' ' is an empty suffix
125 c CALL READ_FLD_XY_RS( bathyFile, ' ', R_low, 0, myThid )
126 C Read the bathymetry using the low-level I/O package
127 c CALL MDSREADFIELD( bathyFile, readBinaryPrec,
128 c & 'RS', 1, R_low, 1, myThid )
129
130 #ifdef ALLOW_MNC
131 ENDIF
132 #endif /* ALLOW_MNC */
133
134
135 ENDIF
136 C- end setup R_low in the interior
137
138 C- fill in the overlap (+ BARRIER):
139 _EXCH_XY_RS(R_low, myThid )
140
141 IF ( debugLevel.GE.debLevB ) THEN
142 c PRINT *, ' Calling plot field', myThid
143 CALL PLOT_FIELD_XYRS( R_low, 'Bottom depths (ini_depths)',
144 & -1, myThid )
145 ENDIF
146 c CALL WRITE_FLD_XY_RS( 'R_low' ,' ', R_low, 0,myThid)
147
148 c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
149
150 C------
151 C 2) Set R_surf = Surface boundary: ocean surface / ground for the atmosphere
152 C------
153
154 IF ( usingPCoords .AND. bathyFile.NE.' ' ) THEN
155 C------ read directly Po_surf from bathyFile (only for backward compatibility)
156
157 CALL READ_REC_XY_RS( bathyFile, Ro_surf, 1, 0, myThid )
158
159 ELSEIF ( topoFile.EQ.' ' ) THEN
160 C------ set default value:
161
162 DO bj = myByLo(myThid), myByHi(myThid)
163 DO bi = myBxLo(myThid), myBxHi(myThid)
164 DO j=1,sNy
165 DO i=1,sNx
166 Ro_surf(i,j,bi,bj) = rF(1)
167 ENDDO
168 ENDDO
169 ENDDO
170 ENDDO
171
172 ELSE
173 C------ read from file:
174
175 C- read surface topography (in m) from topoFile (case topoFile.NE.' '):
176
177 CALL READ_REC_XY_RS( topoFile, topoZ, 1, 0, myThid )
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 This I/O is now done in write_grid.F
191 c CALL WRITE_FLD_XY_RS( 'topo_P',' ',Ro_surf,0,myThid)
192
193 ELSEIF ( buoyancyRelation.EQ.'OCEANICP' ) THEN
194
195 WRITE(msgBuf,'(A,A)') 'S/R INI_DEPTHS: ',
196 & 'from topoFile (in m) to ref.bottom pressure: Not yet coded'
197 CALL PRINT_ERROR( msgBuf , myThid)
198 STOP 'ABNORMAL END: S/R INI_DEPTHS'
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_RS(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) = rF(1)
257 ENDDO
258 ENDDO
259 ENDDO
260 ENDDO
261 ENDIF
262
263 IF ( debugLevel.GE.debLevB ) THEN
264 _BARRIER
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: not necessary since all single-thread IO above are followed
273 C by an EXCH (with BARRIER) + BARRIER within IO routine
274 c _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_RS(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 CcnhSigmaTestingStart
291 CALL SIGMA_TESTING_LOAD_TOPO( myThid )
292 CcnhSigmaTestingEnd
293
294 RETURN
295 END

  ViewVC Help
Powered by ViewVC 1.1.22