/[MITgcm]/MITgcm/verification/advect_xz/code/ini_depths.F
ViewVC logotype

Annotation of /MITgcm/verification/advect_xz/code/ini_depths.F

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


Revision 1.2 - (hide annotations) (download)
Tue Dec 10 03:07:28 2002 UTC (21 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint48f_post, checkpoint51k_post, checkpoint53f_post, checkpoint47j_post, checkpoint54a_pre, checkpoint55c_post, checkpoint53b_pre, checkpoint48d_pre, checkpoint51l_post, checkpoint51j_post, branch-exfmods-tag, checkpoint47e_post, checkpoint47i_post, checkpoint52l_pre, checkpoint48i_post, checkpoint52e_pre, checkpoint57g_pre, checkpoint52j_post, checkpoint47f_post, checkpoint48d_post, checkpoint51o_pre, checkpoint57f_post, checkpoint47c_post, checkpoint50e_post, checkpoint52e_post, checkpoint50c_post, checkpoint51n_pre, checkpoint47d_post, checkpoint57b_post, checkpoint53c_post, checkpoint53d_post, checkpoint57f_pre, checkpoint48a_post, checkpoint55d_pre, checkpoint51f_pre, checkpoint57g_post, checkpoint48e_post, checkpoint57a_post, checkpoint48h_post, checkpoint55j_post, checkpoint56b_post, checkpoint50c_pre, checkpoint57h_pre, branchpoint-genmake2, checkpoint52j_pre, checkpoint54a_post, branch-netcdf, checkpoint50d_pre, checkpoint55h_post, checkpoint51r_post, checkpoint52b_pre, checkpoint52n_post, checkpoint54b_post, checkpoint51i_post, checkpoint57e_post, checkpoint54d_post, checkpoint47h_post, checkpoint48c_post, checkpoint56c_post, checkpoint54e_post, checkpoint55b_post, checkpoint51e_post, checkpoint51b_post, checkpoint51l_pre, checkpoint52m_post, checkpoint51c_post, checkpoint55, checkpoint53a_post, checkpoint55a_post, checkpoint57c_pre, checkpoint48, checkpoint49, checkpoint53b_post, checkpoint55g_post, checkpoint51o_post, checkpoint48g_post, checkpoint57d_post, checkpoint55f_post, checkpoint57i_post, checkpoint51q_post, checkpoint52l_post, checkpoint52k_post, checkpoint57h_post, checkpoint57a_pre, checkpoint54, checkpoint57, checkpoint56, checkpoint51, checkpoint50, checkpoint53, checkpoint52, checkpoint50d_post, checkpoint52d_post, checkpoint51b_pre, checkpoint52a_post, checkpoint57h_done, checkpoint47g_post, checkpoint52b_post, checkpoint53g_post, checkpoint52f_post, checkpoint52c_post, checkpoint51h_pre, checkpoint50g_post, checkpoint51g_post, ecco_c52_e35, checkpoint54f_post, checkpoint51f_post, checkpoint48b_post, checkpoint50b_post, eckpoint57e_pre, checkpoint57c_post, checkpoint50f_post, checkpoint50a_post, checkpoint50f_pre, checkpoint52a_pre, checkpoint47d_pre, checkpoint51d_post, checkpoint48c_pre, checkpoint51m_post, checkpoint51t_post, checkpoint53d_pre, checkpoint55e_post, checkpoint54c_post, checkpoint52i_post, checkpoint51a_post, checkpoint50e_pre, checkpoint50i_post, checkpoint51p_post, checkpoint51n_post, checkpoint55i_post, checkpoint51i_pre, checkpoint52i_pre, checkpoint51u_post, checkpoint52h_pre, checkpoint52f_pre, checkpoint50b_pre, checkpoint56a_post, checkpoint51s_post, checkpoint55d_post
Branch point for: netcdf-sm0, branch-genmake2, branch-nonh, tg2-branch, checkpoint51n_branch, branch-exfmods-curt
Changes since 1.1: +179 -67 lines
 update after changing the standard version

1 jmc 1.2 C $Header: /u/gcmpack/MITgcm/verification/advect_xz/code/ini_depths.F,v 1.1 2001/09/28 02:28:10 adcroft Exp $
2 adcroft 1.1 C $Name: $
3    
4     #include "CPP_OPTIONS.h"
5    
6 jmc 1.2 CBOP
7     C !ROUTINE: INI_DEPTHS
8     C !INTERFACE:
9 adcroft 1.1 SUBROUTINE INI_DEPTHS( myThid )
10 jmc 1.2 C !DESCRIPTION: \bv
11     C *==========================================================*
12     C | SUBROUTINE INI_DEPTHS
13     C | o define R_position of Lower and Surface Boundaries
14     C *==========================================================*
15     C |atmosphere orography:
16     C | define either in term of P_topo or converted from Z_topo
17     C |ocean bathymetry:
18     C | The depths of the bottom of the model is specified in
19     C | terms of an XY map with one depth for each column of
20     C | grid cells. Depths do not have to coincide with the
21     C | model levels. The model lopping algorithm makes it
22     C | possible to represent arbitrary depths.
23     C | The mode depths map also influences the models topology
24     C | By default the model domain wraps around in X and Y.
25     C | This default doubly periodic topology is "supressed"
26     C | if a depth map is defined which closes off all wrap
27     C | around flow.
28     C *==========================================================*
29     C \ev
30    
31     C !USES:
32 adcroft 1.1 IMPLICIT NONE
33     C === Global variables ===
34     #include "SIZE.h"
35     #include "EEPARAMS.h"
36     #include "PARAMS.h"
37     #include "GRID.h"
38 jmc 1.2 #include "SURFACE.h"
39 adcroft 1.1
40 jmc 1.2 C !INPUT/OUTPUT PARAMETERS:
41 adcroft 1.1 C == Routine arguments ==
42     C myThid - Number of this instance of INI_DEPTHS
43     INTEGER myThid
44     CEndOfInterface
45    
46 jmc 1.2 C !LOCAL VARIABLES:
47 adcroft 1.1 C == Local variables ==
48     C iG, jG - Global coordinate index
49     C bi,bj - Loop counters
50     C I,J,K
51     C oldPrec - Temporary used in controlling binary input dataset precision
52 jmc 1.2 C msgBuf - Informational/error meesage buffer
53 adcroft 1.1 INTEGER iG, jG
54     INTEGER bi, bj
55     INTEGER I, J
56 jmc 1.2 CHARACTER*(MAX_LEN_MBUF) msgBuf
57     CEOP
58 adcroft 1.1
59 jmc 1.2 IF (groundAtK1 .AND. bathyFile .NE. ' '
60     & .AND. topoFile .NE. ' ' ) THEN
61     WRITE(msgBuf,'(A,A)')
62     & 'S/R INI_DEPTHS: both bathyFile & topoFile are specified:',
63     & ' select the right one !'
64     CALL PRINT_ERROR( msgBuf , myThid)
65     STOP 'ABNORMAL END: S/R INI_DEPTHS'
66     ENDIF
67    
68     C------
69     C 0) Initialize R_low and Ro_surf (define an empty domain)
70     C------
71     DO bj = myByLo(myThid), myByHi(myThid)
72     DO bi = myBxLo(myThid), myBxHi(myThid)
73     DO j=1-Oly,sNy+Oly
74     DO i=1-Olx,sNx+Olx
75     R_low(i,j,bi,bj) = 0.
76     Ro_surf(i,j,bi,bj) = 0.
77     topoZ(i,j,bi,bj) = 0.
78     ENDDO
79     ENDDO
80     ENDDO
81     ENDDO
82    
83     C------
84     C 1) Set R_low = the Lower (in r sense) boundary of the fluid column :
85     C------
86     IF (groundAtK1 .OR. bathyFile .EQ. ' ') THEN
87     C- e.g., atmosphere : R_low = Top of atmosphere
88     C- ocean : R_low = Bottom
89 adcroft 1.1 DO bj = myByLo(myThid), myByHi(myThid)
90     DO bi = myBxLo(myThid), myBxHi(myThid)
91     DO j=1,sNy
92     DO i=1,sNx
93 jmc 1.2 R_low(i,j,bi,bj) = rF(Nr+1)
94     C-- Specific modif for this experiment (advect_xz):
95     R_low(I,J,bi,bj) = R_low(I,J,bi,bj)
96 adcroft 1.1 & *(1.-0.5*XC(I,J,bi,bj)/(float(Nx)*DelX(1)))
97 jmc 1.2 C-- end of modified part
98 adcroft 1.1 ENDDO
99     ENDDO
100     ENDDO
101     ENDDO
102     ELSE
103 jmc 1.2 _BEGIN_MASTER( myThid )
104 adcroft 1.1 C Read the bathymetry using the mid-level I/O pacakage read_write_rec
105     C The 0 is the "iteration" argument. The 1 is the record number.
106 jmc 1.2 CALL READ_REC_XY_RS( bathyFile, R_low, 1, 0, myThid )
107 adcroft 1.1 C Read the bathymetry using the mid-level I/O pacakage read_write_fld
108     C The 0 is the "iteration" argument. The ' ' is an empty suffix
109 jmc 1.2 c CALL READ_FLD_XY_RS( bathyFile, ' ', R_low, 0, myThid )
110 adcroft 1.1 C Read the bathymetry using the low-level I/O package
111 jmc 1.2 c CALL MDSREADFIELD( bathyFile, readBinaryPrec,
112     c & 'RS', 1, R_low, 1, myThid )
113     _END_MASTER(myThid)
114    
115 adcroft 1.1 ENDIF
116 jmc 1.2 C- end setup R_low in the interior
117 adcroft 1.1
118 jmc 1.2 C- fill in the overlap :
119     _EXCH_XY_R4(R_low, myThid )
120 adcroft 1.1
121 jmc 1.2 c CALL PLOT_FIELD_XYRS(R_low,'Bottom depths (ini_depths)',1,myThid)
122     c _BEGIN_MASTER( myThid )
123     c CALL WRITE_FLD_XY_RS( 'R_low' ,' ', R_low, 0,myThid)
124     c _END_MASTER(myThid)
125    
126     c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
127    
128     C------
129     C 2) Set R_surf = Surface boundary: ocean surface / ground for the atmosphere
130     C------
131    
132     IF ( groundAtK1 .AND. bathyFile.NE.' ' ) THEN
133     C------ read directly Po_surf from bathyFile (only for backward compatibility)
134    
135     _BEGIN_MASTER( myThid )
136     CALL READ_REC_XY_RS( bathyFile, Ro_surf, 1, 0, myThid )
137     _END_MASTER(myThid)
138     _BARRIER
139    
140     ELSEIF ( topoFile.EQ.' ' ) THEN
141     C------ set default value:
142    
143     DO bj = myByLo(myThid), myByHi(myThid)
144     DO bi = myBxLo(myThid), myBxHi(myThid)
145     DO j=1,sNy
146     DO i=1,sNx
147     Ro_surf(i,j,bi,bj) = Ro_SeaLevel
148     ENDDO
149     ENDDO
150     ENDDO
151     ENDDO
152    
153     ELSE
154     C------ read from file:
155    
156     C- read surface topography (in m) from topoFile (case topoFile.NE.' '):
157     _BEGIN_MASTER( myThid )
158     CALL READ_REC_XY_RS( topoFile, topoZ, 1, 0, myThid )
159     _END_MASTER(myThid)
160     _BARRIER
161    
162     IF (buoyancyRelation .EQ. 'ATMOSPHERIC') THEN
163     C----
164     C Convert Surface Geopotential to (reference) Surface Pressure
165     C according to Tref profile, using same discretisation as in calc_phi_hyd
166     C----
167     c _BEGIN_MASTER( myThid )
168     c CALL WRITE_FLD_XY_RS( 'topo_Z',' ',topoZ,0,myThid)
169     c _END_MASTER(myThid)
170    
171     CALL INI_P_GROUND( selectFindRoSurf, topoZ,
172     O Ro_surf,
173     I myThid )
174    
175     _BARRIER
176     _BEGIN_MASTER( myThid )
177     CALL WRITE_FLD_XY_RS( 'topo_P',' ',Ro_surf,0,myThid)
178     _END_MASTER(myThid)
179    
180     ELSE
181     C----
182     C Direct Transfer to Ro_surf :
183     DO bj = myByLo(myThid), myByHi(myThid)
184     DO bi = myBxLo(myThid), myBxHi(myThid)
185     DO j=1,sNy
186     DO i=1,sNx
187     Ro_surf(i,j,bi,bj) = topoZ(i,j,bi,bj)
188     ENDDO
189     ENDDO
190     ENDDO
191     ENDDO
192    
193     ENDIF
194    
195     C------ end case "read topoFile"
196     ENDIF
197    
198     C----- fill in the overlap :
199     _EXCH_XY_R4(Ro_surf, myThid )
200    
201     c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
202    
203     C------
204     C 3) Close the Domain (special configuration).
205     C------
206 adcroft 1.1 IF (groundAtK1) THEN
207     DO bj = myByLo(myThid), myByHi(myThid)
208     DO bi = myBxLo(myThid), myBxHi(myThid)
209     DO j=1-Oly,sNy+Oly
210     DO i=1-Olx,sNx+Olx
211 jmc 1.2 iG = myXGlobalLo-1+(bi-1)*sNx+I
212     jG = myYGlobalLo-1+(bj-1)*sNy+J
213     C Test for eastern edge
214     c IF ( iG .EQ. Nx ) Ro_surf(i,j,bi,bj) = 0.
215     C Test for northern edge
216     c IF ( jG .EQ. Ny ) Ro_surf(i,j,bi,bj) = 0.
217 adcroft 1.1 IF (usingSphericalPolarGrid .AND. abs(yC(I,J,bi,bj)).GE.90. )
218 jmc 1.2 & Ro_surf(I,J,bi,bj) = rF(Nr+1)
219 adcroft 1.1 ENDDO
220     ENDDO
221     ENDDO
222     ENDDO
223     ELSE
224     DO bj = myByLo(myThid), myByHi(myThid)
225     DO bi = myBxLo(myThid), myBxHi(myThid)
226     DO j=1-Oly,sNy+Oly
227     DO i=1-Olx,sNx+Olx
228 jmc 1.2 iG = myXGlobalLo-1+(bi-1)*sNx+I
229     jG = myYGlobalLo-1+(bj-1)*sNy+J
230     C Test for eastern edge
231     c IF ( iG .EQ. Nx ) R_low(i,j,bi,bj) = 0.
232     C Test for northern edge
233     c IF ( jG .EQ. Ny ) R_low(i,j,bi,bj) = 0.
234 adcroft 1.1 IF (usingSphericalPolarGrid .AND. abs(yC(I,J,bi,bj)).GE.90. )
235 jmc 1.2 & R_low(I,J,bi,bj) = Ro_SeaLevel
236 adcroft 1.1 ENDDO
237     ENDDO
238     ENDDO
239     ENDDO
240     ENDIF
241 jmc 1.2
242     c _BEGIN_MASTER( myThid )
243     c CALL WRITE_FLD_XY_RS('Ro_surf',' ',Ro_surf,0,myThid)
244     c _END_MASTER(myThid)
245 adcroft 1.1
246     RETURN
247     END

  ViewVC Help
Powered by ViewVC 1.1.22