/[MITgcm]/MITgcm_contrib/heimbach/OpenAD/code_ad_moc/ini_depths.F
ViewVC logotype

Annotation of /MITgcm_contrib/heimbach/OpenAD/code_ad_moc/ini_depths.F

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


Revision 1.1 - (hide annotations) (download)
Mon Sep 15 22:16:00 2008 UTC (16 years, 10 months ago) by utke
Branch: MAIN
CVS Tags: HEAD
keep this version around

1 utke 1.1 C $Header: /u/gcmpack/MITgcm/model/src/ini_depths.F,v 1.39 2006/11/05 18:36:05 edhill 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     cph IF ( debugLevel.GE.debLevB ) THEN
147     CALL PLOT_FIELD_XYRS( R_low, 'Bottom depths (ini_depths)',
148     & -1, myThid )
149     cph ENDIF
150     c CALL WRITE_FLD_XY_RS( 'R_low' ,' ', R_low, 0,myThid)
151    
152     c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
153    
154     C------
155     C 2) Set R_surf = Surface boundary: ocean surface / ground for the atmosphere
156     C------
157    
158     IF ( usingPCoords .AND. bathyFile.NE.' ' ) THEN
159     C------ read directly Po_surf from bathyFile (only for backward compatibility)
160    
161     CALL READ_REC_XY_RS( bathyFile, Ro_surf, 1, 0, myThid )
162    
163     ELSEIF ( topoFile.EQ.' ' ) THEN
164     C------ set default value:
165    
166     DO bj = myByLo(myThid), myByHi(myThid)
167     DO bi = myBxLo(myThid), myBxHi(myThid)
168     DO j=1,sNy
169     DO i=1,sNx
170     Ro_surf(i,j,bi,bj) = Ro_SeaLevel
171     ENDDO
172     ENDDO
173     ENDDO
174     ENDDO
175    
176     ELSE
177     C------ read from file:
178    
179     C- read surface topography (in m) from topoFile (case topoFile.NE.' '):
180    
181     CALL READ_REC_XY_RS( topoFile, topoZ, 1, 0, myThid )
182     _BARRIER
183    
184     IF (buoyancyRelation .EQ. 'ATMOSPHERIC') THEN
185     C----
186     C Convert Surface Geopotential to (reference) Surface Pressure
187     C according to Tref profile, using same discretisation as in calc_phi_hyd
188     C----
189     c CALL WRITE_FLD_XY_RS( 'topo_Z',' ',topoZ,0,myThid)
190    
191     CALL INI_P_GROUND( 2, topoZ,
192     O Ro_surf,
193     I myThid )
194    
195     c _BARRIER
196     C This I/O is now done in write_grid.F
197     c CALL WRITE_FLD_XY_RS( 'topo_P',' ',Ro_surf,0,myThid)
198    
199     ELSE
200     C----
201     C Direct Transfer to Ro_surf (e.g., to specify upper ocean boundary
202     C below an ice-shelf - NOTE - actually not yet implemented )
203     DO bj = myByLo(myThid), myByHi(myThid)
204     DO bi = myBxLo(myThid), myBxHi(myThid)
205     DO j=1,sNy
206     DO i=1,sNx
207     Ro_surf(i,j,bi,bj) = topoZ(i,j,bi,bj)
208     ENDDO
209     ENDDO
210     ENDDO
211     ENDDO
212    
213     ENDIF
214    
215     C------ end case "read topoFile"
216     ENDIF
217    
218     C----- fill in the overlap (+ BARRIER):
219     _EXCH_XY_R4(Ro_surf, myThid )
220    
221     c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
222    
223     C------
224     C 3) Close the Domain (special configuration).
225     C------
226     IF (usingPCoords) THEN
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 ) Ro_surf(i,j,bi,bj) = 0.
235     C Test for northern edge
236     c IF ( jG .EQ. Ny ) Ro_surf(i,j,bi,bj) = 0.
237     IF (usingSphericalPolarGrid .AND. ABS(yC(I,J,bi,bj)).GE.90. )
238     & Ro_surf(I,J,bi,bj) = rF(Nr+1)
239     ENDDO
240     ENDDO
241     ENDDO
242     ENDDO
243     ELSE
244     DO bj = myByLo(myThid), myByHi(myThid)
245     DO bi = myBxLo(myThid), myBxHi(myThid)
246     DO j=1-Oly,sNy+Oly
247     DO i=1-Olx,sNx+Olx
248     iG = myXGlobalLo-1+(bi-1)*sNx+I
249     jG = myYGlobalLo-1+(bj-1)*sNy+J
250     C Test for eastern edge
251     c IF ( iG .EQ. Nx ) R_low(i,j,bi,bj) = 0.
252     C Test for northern edge
253     c IF ( jG .EQ. Ny ) R_low(i,j,bi,bj) = 0.
254     IF (usingSphericalPolarGrid .AND. ABS(yC(I,J,bi,bj)).GE.90. )
255     & R_low(I,J,bi,bj) = Ro_SeaLevel
256     ENDDO
257     ENDDO
258     ENDDO
259     ENDDO
260     ENDIF
261    
262     IF ( debugLevel.GE.debLevB ) THEN
263     CALL PLOT_FIELD_XYRS( Ro_surf,
264     & 'Surface reference r-position (ini_depths)',
265     & -1, myThid )
266     ENDIF
267     c CALL WRITE_FLD_XY_RS('Ro_surf',' ',Ro_surf,0,myThid)
268    
269     #ifdef ALLOW_SHELFICE
270     c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
271     IF ( useShelfIce ) THEN
272     C------
273     C 4) Set R_shelfIce = the Lower (in r sense) boundary of floating shelfice :
274     C------
275     IF (usingPCoords .OR. shelfIceFile .EQ. ' ') THEN
276     C- e.g., atmosphere : R_low = Top of atmosphere
277     C- ocean : R_low = Bottom
278     DO bj = myByLo(myThid), myByHi(myThid)
279     DO bi = myBxLo(myThid), myBxHi(myThid)
280     DO j=1,sNy
281     DO i=1,sNx
282     R_shelfIce(i,j,bi,bj) = 0. _d 0
283     ENDDO
284     ENDDO
285     ENDDO
286     ENDDO
287     ELSE
288     C Read the shelfIce draught using the mid-level I/O pacakage read_write_rec
289     C The 0 is the "iteration" argument. The 1 is the record number.
290     CALL READ_REC_XY_RS( shelfIceFile, R_shelfIce, 1, 0, myThid )
291     C Read the shelfIce draught using the mid-level I/O pacakage read_write_fld
292     C The 0 is the "iteration" argument. The ' ' is an empty suffix
293     C CALL READ_FLD_XY_RS( shelfIceFile, ' ', R_shelfIce, 0, myThid )
294     c CALL READ_FLD_XY_RS( bathyFile, ' ', R_low, 0, myThid )
295     C Read the selfIce draught using the low-level I/O package
296     c CALL MDSREADFIELD( shelfIceFile, readBinaryPrec,
297     c & 'RS', 1, R_selfIce, 1, myThid )
298    
299     ENDIF
300     C- end setup R_shelfIce in the interior
301    
302     C- fill in the overlap (+ BARRIER):
303     _EXCH_XY_R4(R_shelfIce, myThid )
304    
305     c CALL PLOT_FIELD_XYRS(R_selfIce,'Shelf ice draught (ini_depths)',
306     c & 1,myThid)
307     CML CALL WRITE_FLD_XY_RS( 'R_shelfIce' ,' ', R_shelfIce, 0,myThid)
308     ENDIF
309     #endif /* ALLOW_SHELFICE */
310    
311     C-- Everyone else must wait for the depth to be loaded
312     C- note: might not be necessary since all single-thread IO above
313     C are followed by an EXCH (with BARRIER)
314     _BARRIER
315    
316     RETURN
317     END

  ViewVC Help
Powered by ViewVC 1.1.22