/[MITgcm]/MITgcm/verification/aim.5l_Equatorial_Channel/code/ini_depths.F
ViewVC logotype

Contents of /MITgcm/verification/aim.5l_Equatorial_Channel/code/ini_depths.F

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


Revision 1.12 - (show annotations) (download)
Tue Jan 27 15:38:07 2009 UTC (15 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61l, checkpoint61j, checkpoint61k, checkpoint61i
Changes since 1.11: +3 -3 lines
rename thetaMin,phiMin -> xgOrigin,ygOrigin

1 C $Header: /u/gcmpack/MITgcm/verification/aim.5l_Equatorial_Channel/code/ini_depths.F,v 1.11 2006/10/18 20:44:37 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
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 #ifdef ALLOW_SHELFICE
83 R_shelfIce(i,j,bi,bj) = 0. _d 0
84 #endif /* ALLOW_SHELFICE */
85 ENDDO
86 ENDDO
87 ENDDO
88 ENDDO
89
90 C- Need to synchronize here before doing master-thread IO
91 _BARRIER
92
93 C------
94 C 1) Set R_low = the Lower (in r sense) boundary of the fluid column :
95 C------
96 IF (usingPCoords .OR. bathyFile .EQ. ' ') THEN
97 C- e.g., atmosphere : R_low = Top of atmosphere
98 C- ocean : R_low = Bottom
99 DO bj = myByLo(myThid), myByHi(myThid)
100 DO bi = myBxLo(myThid), myBxHi(myThid)
101 DO j=1,sNy
102 DO i=1,sNx
103 R_low(i,j,bi,bj) = rF(Nr+1)
104 ENDDO
105 ENDDO
106 ENDDO
107 ENDDO
108 ELSE
109 C Read the bathymetry using the mid-level I/O package read_write_rec
110 C The 0 is the "iteration" argument. The 1 is the record number.
111 CALL READ_REC_XY_RS( bathyFile, R_low, 1, 0, myThid )
112 C Read the bathymetry using the mid-level I/O package read_write_fld
113 C The 0 is the "iteration" argument. The ' ' is an empty suffix
114 c CALL READ_FLD_XY_RS( bathyFile, ' ', R_low, 0, myThid )
115 C Read the bathymetry using the low-level I/O package
116 c CALL MDSREADFIELD( bathyFile, readBinaryPrec,
117 c & 'RS', 1, R_low, 1, myThid )
118 ENDIF
119 C- end setup R_low in the interior
120
121 C- fill in the overlap (+ BARRIER):
122 _EXCH_XY_R4(R_low, myThid )
123
124 IF ( debugLevel.GE.debLevB ) THEN
125 c PRINT *, ' Calling plot field', myThid
126 CALL PLOT_FIELD_XYRS( R_low, 'Bottom depths (ini_depths)',
127 & -1, myThid )
128 ENDIF
129 c CALL WRITE_FLD_XY_RS( 'R_low' ,' ', R_low, 0,myThid)
130
131 c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
132
133 C------
134 C 2) Set R_surf = Surface boundary: ocean surface / ground for the atmosphere
135 C------
136
137 IF ( usingPCoords .AND. bathyFile.NE.' ' ) THEN
138 C------ read directly Po_surf from bathyFile (only for backward compatibility)
139
140 CALL READ_REC_XY_RS( bathyFile, Ro_surf, 1, 0, myThid )
141
142 ELSEIF ( topoFile.EQ.' ' ) THEN
143 C------ set default value:
144
145 DO bj = myByLo(myThid), myByHi(myThid)
146 DO bi = myBxLo(myThid), myBxHi(myThid)
147 DO j=1,sNy
148 DO i=1,sNx
149 Ro_surf(i,j,bi,bj) = Ro_SeaLevel
150 ENDDO
151 ENDDO
152 ENDDO
153 ENDDO
154
155 ELSE
156 C------ read from file:
157
158 C- read surface topography (in m) from topoFile (case topoFile.NE.' '):
159
160 CALL READ_REC_XY_RS( topoFile, topoZ, 1, 0, myThid )
161 _BARRIER
162
163 IF (buoyancyRelation .EQ. 'ATMOSPHERIC') THEN
164 C----
165 C Convert Surface Geopotential to (reference) Surface Pressure
166 C according to Tref profile, using same discretisation as in calc_phi_hyd
167 C----
168 c CALL WRITE_FLD_XY_RS( 'topo_Z',' ',topoZ,0,myThid)
169
170 CALL INI_P_GROUND( 2, topoZ,
171 O Ro_surf,
172 I myThid )
173
174 c _BARRIER
175 C This I/O is now done in write_grid.F
176 c CALL WRITE_FLD_XY_RS( 'topo_P',' ',Ro_surf,0,myThid)
177
178 ELSE
179 C----
180 C Direct Transfer to Ro_surf (e.g., to specify upper ocean boundary
181 C below an ice-shelf - NOTE - actually not yet implemented )
182 DO bj = myByLo(myThid), myByHi(myThid)
183 DO bi = myBxLo(myThid), myBxHi(myThid)
184 DO j=1,sNy
185 DO i=1,sNx
186 Ro_surf(i,j,bi,bj) = topoZ(i,j,bi,bj)
187 ENDDO
188 ENDDO
189 ENDDO
190 ENDDO
191
192 ENDIF
193
194 C------ end case "read topoFile"
195 ENDIF
196
197 C----- fill in the overlap (+ BARRIER):
198 _EXCH_XY_R4(Ro_surf, myThid )
199
200 c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
201
202 C------
203 C 3) Close the Domain (special configuration).
204 C------
205 IF (usingPCoords) THEN
206 DO bj = myByLo(myThid), myByHi(myThid)
207 DO bi = myBxLo(myThid), myBxHi(myThid)
208 DO j=1-Oly,sNy+Oly
209 DO i=1-Olx,sNx+Olx
210 iG = myXGlobalLo-1+(bi-1)*sNx+I
211 jG = myYGlobalLo-1+(bj-1)*sNy+J
212 C Test for eastern edge
213 c IF ( iG .EQ. Nx ) Ro_surf(i,j,bi,bj) = 0.
214 C Test for northern edge
215 c IF ( jG .EQ. Ny ) Ro_surf(i,j,bi,bj) = 0.
216 c- Domain : Symetric % Eq. & closed at N & S boundaries:
217 IF (usingSphericalPolarGrid .AND.
218 & ABS(yC(I,J,bi,bj)).GE.-ygOrigin)
219 & Ro_surf(I,J,bi,bj) = rF(Nr+1)
220 IF (usingSphericalPolarGrid .AND. ABS(yC(I,J,bi,bj)).GE.90. )
221 & Ro_surf(I,J,bi,bj) = rF(Nr+1)
222 ENDDO
223 ENDDO
224 ENDDO
225 ENDDO
226 ELSE
227 DO bj = myByLo(myThid), myByHi(myThid)
228 DO bi = myBxLo(myThid), myBxHi(myThid)
229 DO j=1-Oly,sNy+Oly
230 DO i=1-Olx,sNx+Olx
231 iG = myXGlobalLo-1+(bi-1)*sNx+I
232 jG = myYGlobalLo-1+(bj-1)*sNy+J
233 C Test for eastern edge
234 c IF ( iG .EQ. Nx ) R_low(i,j,bi,bj) = 0.
235 C Test for northern edge
236 c IF ( jG .EQ. Ny ) R_low(i,j,bi,bj) = 0.
237 c- Domain : Symetric % Eq. & closed at N & S boundaries:
238 IF (usingSphericalPolarGrid .AND.
239 & ABS(yC(I,J,bi,bj)).GE.-ygOrigin)
240 & R_low(I,J,bi,bj) = Ro_SeaLevel
241 IF (usingSphericalPolarGrid .AND. ABS(yC(I,J,bi,bj)).GE.90. )
242 & R_low(I,J,bi,bj) = Ro_SeaLevel
243 ENDDO
244 ENDDO
245 ENDDO
246 ENDDO
247 ENDIF
248
249 IF ( debugLevel.GE.debLevB ) THEN
250 CALL PLOT_FIELD_XYRS( Ro_surf,
251 & 'Surface reference r-position (ini_depths)',
252 & -1, myThid )
253 ENDIF
254 c CALL WRITE_FLD_XY_RS('Ro_surf',' ',Ro_surf,0,myThid)
255
256 #ifdef ALLOW_SHELFICE
257 c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
258 IF ( useShelfIce ) THEN
259 C------
260 C 4) Set R_shelfIce = the Lower (in r sense) boundary of floating shelfice :
261 C------
262 IF (usingPCoords .OR. shelfIceFile .EQ. ' ') THEN
263 C- e.g., atmosphere : R_low = Top of atmosphere
264 C- ocean : R_low = Bottom
265 DO bj = myByLo(myThid), myByHi(myThid)
266 DO bi = myBxLo(myThid), myBxHi(myThid)
267 DO j=1,sNy
268 DO i=1,sNx
269 R_shelfIce(i,j,bi,bj) = 0. _d 0
270 ENDDO
271 ENDDO
272 ENDDO
273 ENDDO
274 ELSE
275 C Read the shelfIce draught using the mid-level I/O pacakage read_write_rec
276 C The 0 is the "iteration" argument. The 1 is the record number.
277 CALL READ_REC_XY_RS( shelfIceFile, R_shelfIce, 1, 0, myThid )
278 C Read the shelfIce draught using the mid-level I/O pacakage read_write_fld
279 C The 0 is the "iteration" argument. The ' ' is an empty suffix
280 C CALL READ_FLD_XY_RS( shelfIceFile, ' ', R_shelfIce, 0, myThid )
281 c CALL READ_FLD_XY_RS( bathyFile, ' ', R_low, 0, myThid )
282 C Read the selfIce draught using the low-level I/O package
283 c CALL MDSREADFIELD( shelfIceFile, readBinaryPrec,
284 c & 'RS', 1, R_selfIce, 1, myThid )
285
286 ENDIF
287 C- end setup R_shelfIce in the interior
288
289 C- fill in the overlap (+ BARRIER):
290 _EXCH_XY_R4(R_shelfIce, myThid )
291
292 c CALL PLOT_FIELD_XYRS(R_selfIce,'Shelf ice draught (ini_depths)',
293 c & 1,myThid)
294 CML CALL WRITE_FLD_XY_RS( 'R_shelfIce' ,' ', R_shelfIce, 0,myThid)
295 ENDIF
296 #endif /* ALLOW_SHELFICE */
297
298 C-- Everyone else must wait for the depth to be loaded
299 C- note: might not be necessary since all single-thread IO above
300 C are followed by an EXCH (with BARRIER)
301 _BARRIER
302
303 RETURN
304 END

  ViewVC Help
Powered by ViewVC 1.1.22