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

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

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


Revision 1.37 - (hide annotations) (download)
Sat Sep 2 23:36:16 2006 UTC (17 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58p_post
Changes since 1.36: +12 -4 lines
call PLOT_FIELD R_low & Ro_surf only if debugLevel >= 2 (since there is
 a 2nd call in ini_mask_etc)

1 jmc 1.37 C $Header: /u/gcmpack/MITgcm/model/src/ini_depths.F,v 1.36 2006/08/04 01:18:48 jmc Exp $
2 adcroft 1.23 C $Name: $
3 cnh 1.1
4 mlosch 1.35 #include "PACKAGES_CONFIG.h"
5 cnh 1.16 #include "CPP_OPTIONS.h"
6 cnh 1.1
7 cnh 1.26 CBOP
8     C !ROUTINE: INI_DEPTHS
9     C !INTERFACE:
10 cnh 1.1 SUBROUTINE INI_DEPTHS( myThid )
11 cnh 1.26 C !DESCRIPTION: \bv
12     C *==========================================================*
13 jmc 1.31 C | SUBROUTINE INI_DEPTHS
14     C | o define R_position of Lower and Surface Boundaries
15 cnh 1.26 C *==========================================================*
16 jmc 1.31 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 cnh 1.26 C *==========================================================*
30     C \ev
31    
32     C !USES:
33 adcroft 1.17 IMPLICIT NONE
34 cnh 1.1 C === Global variables ===
35     #include "SIZE.h"
36     #include "EEPARAMS.h"
37     #include "PARAMS.h"
38     #include "GRID.h"
39 jmc 1.27 #include "SURFACE.h"
40 mlosch 1.35 #ifdef ALLOW_SHELFICE
41     # include "SHELFICE.h"
42     #endif /* ALLOW_SHELFICE */
43 cnh 1.1
44 cnh 1.26 C !INPUT/OUTPUT PARAMETERS:
45 cnh 1.1 C == Routine arguments ==
46     C myThid - Number of this instance of INI_DEPTHS
47     INTEGER myThid
48     CEndOfInterface
49    
50 cnh 1.26 C !LOCAL VARIABLES:
51 cnh 1.1 C == Local variables ==
52     C iG, jG - Global coordinate index
53 jmc 1.29 C bi,bj - Tile indices
54     C I,J,K - Loop counters
55 cnh 1.13 C oldPrec - Temporary used in controlling binary input dataset precision
56 jmc 1.31 C msgBuf - Informational/error meesage buffer
57 cnh 1.1 INTEGER iG, jG
58     INTEGER bi, bj
59 adcroft 1.23 INTEGER I, J
60 jmc 1.24 CHARACTER*(MAX_LEN_MBUF) msgBuf
61 cnh 1.26 CEOP
62 jmc 1.24
63 jmc 1.31 IF (usingPCoords .AND. bathyFile .NE. ' '
64     & .AND. topoFile .NE. ' ' ) THEN
65 jmc 1.24 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 jmc 1.36 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 jmc 1.24 ENDDO
86     ENDDO
87     ENDDO
88     ENDDO
89 cnh 1.1
90 jmc 1.36 C- Need to synchronize here before doing master-thread IO
91     _BARRIER
92    
93 jmc 1.24 C------
94     C 1) Set R_low = the Lower (in r sense) boundary of the fluid column :
95     C------
96 jmc 1.31 IF (usingPCoords .OR. bathyFile .EQ. ' ') THEN
97 jmc 1.24 C- e.g., atmosphere : R_low = Top of atmosphere
98 jmc 1.25 C- ocean : R_low = Bottom
99 cnh 1.4 DO bj = myByLo(myThid), myByHi(myThid)
100     DO bi = myBxLo(myThid), myBxHi(myThid)
101 jmc 1.25 DO j=1,sNy
102     DO i=1,sNx
103 jmc 1.24 R_low(i,j,bi,bj) = rF(Nr+1)
104 cnh 1.4 ENDDO
105 cnh 1.1 ENDDO
106     ENDDO
107     ENDDO
108 cnh 1.4 ELSE
109 jmc 1.36 C Read the bathymetry using the mid-level I/O package read_write_rec
110 adcroft 1.19 C The 0 is the "iteration" argument. The 1 is the record number.
111 jmc 1.24 CALL READ_REC_XY_RS( bathyFile, R_low, 1, 0, myThid )
112 jmc 1.36 C Read the bathymetry using the mid-level I/O package read_write_fld
113 adcroft 1.19 C The 0 is the "iteration" argument. The ' ' is an empty suffix
114 jmc 1.24 c CALL READ_FLD_XY_RS( bathyFile, ' ', R_low, 0, myThid )
115 adcroft 1.19 C Read the bathymetry using the low-level I/O package
116 jmc 1.24 c CALL MDSREADFIELD( bathyFile, readBinaryPrec,
117     c & 'RS', 1, R_low, 1, myThid )
118 jmc 1.25 ENDIF
119     C- end setup R_low in the interior
120 jmc 1.24
121 jmc 1.36 C- fill in the overlap (+ BARRIER):
122 jmc 1.25 _EXCH_XY_R4(R_low, myThid )
123 adcroft 1.21
124 jmc 1.37 IF ( debugLevel.GE.debLevB ) THEN
125 heimbach 1.34 c PRINT *, ' Calling plot field', myThid
126 jmc 1.37 CALL PLOT_FIELD_XYRS( R_low, 'Bottom depths (ini_depths)',
127     & -1, myThid )
128     ENDIF
129 jmc 1.31 c CALL WRITE_FLD_XY_RS( 'R_low' ,' ', R_low, 0,myThid)
130 jmc 1.24
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 jmc 1.31 IF ( usingPCoords .AND. bathyFile.NE.' ' ) THEN
138     C------ read directly Po_surf from bathyFile (only for backward compatibility)
139 jmc 1.25
140     CALL READ_REC_XY_RS( bathyFile, Ro_surf, 1, 0, myThid )
141    
142     ELSEIF ( topoFile.EQ.' ' ) THEN
143     C------ set default value:
144 jmc 1.24
145     DO bj = myByLo(myThid), myByHi(myThid)
146     DO bi = myBxLo(myThid), myBxHi(myThid)
147 jmc 1.25 DO j=1,sNy
148     DO i=1,sNx
149 jmc 1.24 Ro_surf(i,j,bi,bj) = Ro_SeaLevel
150     ENDDO
151     ENDDO
152     ENDDO
153     ENDDO
154    
155     ELSE
156 jmc 1.25 C------ read from file:
157 jmc 1.24
158     C- read surface topography (in m) from topoFile (case topoFile.NE.' '):
159 jmc 1.36
160 jmc 1.27 CALL READ_REC_XY_RS( topoFile, topoZ, 1, 0, myThid )
161 jmc 1.24 _BARRIER
162    
163     IF (buoyancyRelation .EQ. 'ATMOSPHERIC') THEN
164 jmc 1.25 C----
165 jmc 1.31 C Convert Surface Geopotential to (reference) Surface Pressure
166 jmc 1.24 C according to Tref profile, using same discretisation as in calc_phi_hyd
167 jmc 1.25 C----
168 jmc 1.31 c CALL WRITE_FLD_XY_RS( 'topo_Z',' ',topoZ,0,myThid)
169 jmc 1.24
170 jmc 1.29 CALL INI_P_GROUND( 2, topoZ,
171     O Ro_surf,
172 jmc 1.28 I myThid )
173 jmc 1.24
174 jmc 1.36 c _BARRIER
175 adcroft 1.30 C This I/O is now done in write_grid.F
176 jmc 1.31 c CALL WRITE_FLD_XY_RS( 'topo_P',' ',Ro_surf,0,myThid)
177 jmc 1.24
178     ELSE
179 jmc 1.25 C----
180 jmc 1.31 C Direct Transfer to Ro_surf (e.g., to specify upper ocean boundary
181 jmc 1.29 C below an ice-shelf - NOTE - actually not yet implemented )
182 jmc 1.24 DO bj = myByLo(myThid), myByHi(myThid)
183     DO bi = myBxLo(myThid), myBxHi(myThid)
184     DO j=1,sNy
185     DO i=1,sNx
186 jmc 1.27 Ro_surf(i,j,bi,bj) = topoZ(i,j,bi,bj)
187 jmc 1.24 ENDDO
188     ENDDO
189     ENDDO
190     ENDDO
191 jmc 1.25
192 jmc 1.24 ENDIF
193    
194 jmc 1.25 C------ end case "read topoFile"
195     ENDIF
196 jmc 1.24
197 jmc 1.36 C----- fill in the overlap (+ BARRIER):
198 jmc 1.25 _EXCH_XY_R4(Ro_surf, myThid )
199 jmc 1.24
200     c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
201 adcroft 1.21
202 jmc 1.24 C------
203     C 3) Close the Domain (special configuration).
204     C------
205 jmc 1.31 IF (usingPCoords) THEN
206 adcroft 1.21 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 jmc 1.24 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 jmc 1.37 IF (usingSphericalPolarGrid .AND. ABS(yC(I,J,bi,bj)).GE.90. )
217 jmc 1.24 & Ro_surf(I,J,bi,bj) = rF(Nr+1)
218 adcroft 1.23 ENDDO
219     ENDDO
220     ENDDO
221     ENDDO
222     ELSE
223     DO bj = myByLo(myThid), myByHi(myThid)
224     DO bi = myBxLo(myThid), myBxHi(myThid)
225     DO j=1-Oly,sNy+Oly
226     DO i=1-Olx,sNx+Olx
227 jmc 1.24 iG = myXGlobalLo-1+(bi-1)*sNx+I
228     jG = myYGlobalLo-1+(bj-1)*sNy+J
229     C Test for eastern edge
230     c IF ( iG .EQ. Nx ) R_low(i,j,bi,bj) = 0.
231     C Test for northern edge
232     c IF ( jG .EQ. Ny ) R_low(i,j,bi,bj) = 0.
233 jmc 1.37 IF (usingSphericalPolarGrid .AND. ABS(yC(I,J,bi,bj)).GE.90. )
234 jmc 1.24 & R_low(I,J,bi,bj) = Ro_SeaLevel
235 adcroft 1.21 ENDDO
236     ENDDO
237     ENDDO
238     ENDDO
239     ENDIF
240 jmc 1.24
241 jmc 1.37 IF ( debugLevel.GE.debLevB ) THEN
242     CALL PLOT_FIELD_XYRS( Ro_surf,
243     & 'Surface reference r-position (ini_depths)',
244     & -1, myThid )
245     ENDIF
246 jmc 1.31 c CALL WRITE_FLD_XY_RS('Ro_surf',' ',Ro_surf,0,myThid)
247 adcroft 1.23
248 mlosch 1.35 #ifdef ALLOW_SHELFICE
249     c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
250     IF ( useShelfIce ) THEN
251     C------
252     C 4) Set R_shelfIce = the Lower (in r sense) boundary of floating shelfice :
253     C------
254     IF (usingPCoords .OR. shelfIceFile .EQ. ' ') THEN
255     C- e.g., atmosphere : R_low = Top of atmosphere
256     C- ocean : R_low = Bottom
257     DO bj = myByLo(myThid), myByHi(myThid)
258     DO bi = myBxLo(myThid), myBxHi(myThid)
259     DO j=1,sNy
260     DO i=1,sNx
261     R_shelfIce(i,j,bi,bj) = 0. _d 0
262     ENDDO
263     ENDDO
264     ENDDO
265     ENDDO
266     ELSE
267     C Read the shelfIce draught using the mid-level I/O pacakage read_write_rec
268     C The 0 is the "iteration" argument. The 1 is the record number.
269     CALL READ_REC_XY_RS( shelfIceFile, R_shelfIce, 1, 0, myThid )
270     C Read the shelfIce draught using the mid-level I/O pacakage read_write_fld
271     C The 0 is the "iteration" argument. The ' ' is an empty suffix
272     C CALL READ_FLD_XY_RS( shelfIceFile, ' ', R_shelfIce, 0, myThid )
273     c CALL READ_FLD_XY_RS( bathyFile, ' ', R_low, 0, myThid )
274     C Read the selfIce draught using the low-level I/O package
275     c CALL MDSREADFIELD( shelfIceFile, readBinaryPrec,
276     c & 'RS', 1, R_selfIce, 1, myThid )
277    
278     ENDIF
279     C- end setup R_shelfIce in the interior
280    
281 jmc 1.36 C- fill in the overlap (+ BARRIER):
282 mlosch 1.35 _EXCH_XY_R4(R_shelfIce, myThid )
283    
284     c CALL PLOT_FIELD_XYRS(R_selfIce,'Shelf ice draught (ini_depths)',
285     c & 1,myThid)
286     CML CALL WRITE_FLD_XY_RS( 'R_shelfIce' ,' ', R_shelfIce, 0,myThid)
287     ENDIF
288     #endif /* ALLOW_SHELFICE */
289    
290 cnh 1.1 RETURN
291     END

  ViewVC Help
Powered by ViewVC 1.1.22