/[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.26 - (hide annotations) (download)
Wed Sep 26 18:09:15 2001 UTC (22 years, 8 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint46n_post, checkpoint44e_post, checkpoint46l_post, checkpoint46g_pre, release1_p13_pre, checkpoint46f_post, checkpoint44f_post, checkpoint46b_post, checkpoint43a-release1mods, release1_p13, checkpoint46l_pre, chkpt44d_post, release1_p8, release1_p9, release1_p1, release1_p2, release1_p3, release1_p4, release1_p5, release1_p6, release1_p7, checkpoint44e_pre, release1_b1, checkpoint43, release1_chkpt44d_post, release1_p11, icebear5, icebear4, icebear3, icebear2, checkpoint46d_pre, release1-branch_tutorials, checkpoint45d_post, checkpoint46j_pre, chkpt44a_post, checkpoint44h_pre, checkpoint46a_post, checkpoint46j_post, checkpoint46k_post, chkpt44c_pre, checkpoint45a_post, ecco_c44_e19, ecco_c44_e18, ecco_c44_e17, ecco_c44_e16, release1_p12, release1_p10, release1_p16, release1_p17, release1_p14, release1_p15, checkpoint44g_post, checkpoint46e_pre, checkpoint45b_post, checkpoint46b_pre, release1-branch-end, release1_final_v1, checkpoint46c_pre, checkpoint46, checkpoint44b_post, checkpoint46h_pre, checkpoint46m_post, checkpoint46a_pre, checkpoint45c_post, ecco_ice2, ecco_ice1, checkpoint44h_post, checkpoint46g_post, release1_p12_pre, ecco_c44_e22, ecco_c44_e25, chkpt44a_pre, checkpoint46i_post, ecco_c44_e23, ecco_c44_e20, ecco_c44_e21, ecco_c44_e26, ecco_c44_e27, ecco_c44_e24, checkpoint46c_post, ecco-branch-mod1, ecco-branch-mod2, ecco-branch-mod3, ecco-branch-mod4, ecco-branch-mod5, checkpoint46e_post, release1_beta1, checkpoint44b_pre, checkpoint42, checkpoint41, checkpoint47, checkpoint44, checkpoint45, checkpoint46h_post, chkpt44c_post, checkpoint44f_pre, checkpoint46d_post, release1-branch_branchpoint
Branch point for: c24_e25_ice, release1_final, release1-branch, release1, ecco-branch, release1_50yr, icebear, release1_coupled
Changes since 1.25: +29 -22 lines
Bringing comments up to data and formatting for document extraction.

1 cnh 1.26 C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_depths.F,v 1.25 2001/07/09 16:46:29 jmc Exp $
2 adcroft 1.23 C $Name: $
3 cnh 1.1
4 cnh 1.16 #include "CPP_OPTIONS.h"
5 cnh 1.1
6 cnh 1.26 CBOP
7     C !ROUTINE: INI_DEPTHS
8     C !INTERFACE:
9 cnh 1.1 SUBROUTINE INI_DEPTHS( myThid )
10 cnh 1.26 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.17 IMPLICIT NONE
33 cnh 1.1 C === Global variables ===
34     #include "SIZE.h"
35     #include "EEPARAMS.h"
36     #include "PARAMS.h"
37     #include "GRID.h"
38    
39 cnh 1.26 C !INPUT/OUTPUT PARAMETERS:
40 cnh 1.1 C == Routine arguments ==
41     C myThid - Number of this instance of INI_DEPTHS
42     INTEGER myThid
43     CEndOfInterface
44    
45 cnh 1.26 C !LOCAL VARIABLES:
46 jmc 1.24 C == Local variables in common ==
47     C Hloc - Temporary array used to read surface topography
48     C has to be in common for multi threading
49     COMMON / LOCAL_INI_DEPTHS / Hloc
50     _RS Hloc(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
51 cnh 1.1 C == Local variables ==
52     C iG, jG - Global coordinate index
53     C bi,bj - Loop counters
54 adcroft 1.21 C I,J,K
55 cnh 1.13 C oldPrec - Temporary used in controlling binary input dataset precision
56 jmc 1.24 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     IF (groundAtK1 .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.
80     Ro_surf(i,j,bi,bj) = 0.
81     ENDDO
82     ENDDO
83     ENDDO
84     ENDDO
85 cnh 1.1
86 jmc 1.24 C------
87     C 1) Set R_low = the Lower (in r sense) boundary of the fluid column :
88     C------
89 jmc 1.25 IF (groundAtK1 .OR. bathyFile .EQ. ' ') THEN
90 jmc 1.24 C- e.g., atmosphere : R_low = Top of atmosphere
91 jmc 1.25 C- ocean : R_low = Bottom
92 cnh 1.4 DO bj = myByLo(myThid), myByHi(myThid)
93     DO bi = myBxLo(myThid), myBxHi(myThid)
94 jmc 1.25 DO j=1,sNy
95     DO i=1,sNx
96 jmc 1.24 R_low(i,j,bi,bj) = rF(Nr+1)
97 cnh 1.4 ENDDO
98 cnh 1.1 ENDDO
99     ENDDO
100     ENDDO
101 cnh 1.4 ELSE
102 jmc 1.24 _BEGIN_MASTER( myThid )
103 adcroft 1.19 C Read the bathymetry using the mid-level I/O pacakage read_write_rec
104     C The 0 is the "iteration" argument. The 1 is the record number.
105 jmc 1.24 CALL READ_REC_XY_RS( bathyFile, R_low, 1, 0, myThid )
106 adcroft 1.19 C Read the bathymetry using the mid-level I/O pacakage read_write_fld
107     C The 0 is the "iteration" argument. The ' ' is an empty suffix
108 jmc 1.24 c CALL READ_FLD_XY_RS( bathyFile, ' ', R_low, 0, myThid )
109 adcroft 1.19 C Read the bathymetry using the low-level I/O package
110 jmc 1.24 c CALL MDSREADFIELD( bathyFile, readBinaryPrec,
111     c & 'RS', 1, R_low, 1, myThid )
112     _END_MASTER(myThid)
113    
114 jmc 1.25 ENDIF
115     C- end setup R_low in the interior
116 jmc 1.24
117 jmc 1.25 C- fill in the overlap :
118     _EXCH_XY_R4(R_low, myThid )
119 adcroft 1.21
120 jmc 1.24 c CALL PLOT_FIELD_XYRS(R_low,'Bottom depths (ini_depths)',1,myThid)
121     c _BEGIN_MASTER( myThid )
122     c CALL WRITE_FLD_XY_RS( 'R_low' ,' ', R_low, 0,myThid)
123     c _END_MASTER(myThid)
124    
125     c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
126    
127     C------
128     C 2) Set R_surf = Surface boundary: ocean surface / ground for the atmosphere
129     C------
130    
131 jmc 1.25 IF ( groundAtK1 .AND. bathyFile.NE.' ' ) THEN
132     C------ read directly Po_surf from bathyFile (only for backward compatibility)
133    
134     _BEGIN_MASTER( myThid )
135     CALL READ_REC_XY_RS( bathyFile, Ro_surf, 1, 0, myThid )
136     _END_MASTER(myThid)
137     _BARRIER
138    
139     ELSEIF ( topoFile.EQ.' ' ) THEN
140     C------ set default value:
141 jmc 1.24
142     DO bj = myByLo(myThid), myByHi(myThid)
143     DO bi = myBxLo(myThid), myBxHi(myThid)
144 jmc 1.25 DO j=1,sNy
145     DO i=1,sNx
146 jmc 1.24 Ro_surf(i,j,bi,bj) = Ro_SeaLevel
147     ENDDO
148     ENDDO
149     ENDDO
150     ENDDO
151    
152     ELSE
153 jmc 1.25 C------ read from file:
154 jmc 1.24
155     C- read surface topography (in m) from topoFile (case topoFile.NE.' '):
156     _BEGIN_MASTER( myThid )
157     CALL READ_REC_XY_RS( topoFile, Hloc, 1, 0, myThid )
158     _END_MASTER(myThid)
159     _BARRIER
160    
161     IF (buoyancyRelation .EQ. 'ATMOSPHERIC') THEN
162 jmc 1.25 C----
163 jmc 1.24 C Convert Surface Geopotential to (reference) Surface Pressure
164     C according to Tref profile, using same discretisation as in calc_phi_hyd
165 jmc 1.25 C----
166 jmc 1.24 c _BEGIN_MASTER( myThid )
167     c CALL WRITE_FLD_XY_RS( 'topo_Z',' ',Hloc,0,myThid)
168     c _END_MASTER(myThid)
169    
170     CALL INI_P_GROUND( Hloc, Ro_surf, myThid )
171    
172     _BARRIER
173     _BEGIN_MASTER( myThid )
174     CALL WRITE_FLD_XY_RS( 'topo_P',' ',Ro_surf,0,myThid)
175     _END_MASTER(myThid)
176    
177     ELSE
178 jmc 1.25 C----
179 jmc 1.24 C Direct Transfer to Ro_surf :
180     DO bj = myByLo(myThid), myByHi(myThid)
181     DO bi = myBxLo(myThid), myBxHi(myThid)
182     DO j=1,sNy
183     DO i=1,sNx
184     Ro_surf(i,j,bi,bj) = Hloc(i,j,bi,bj)
185     ENDDO
186     ENDDO
187     ENDDO
188     ENDDO
189 jmc 1.25
190 jmc 1.24 ENDIF
191    
192 jmc 1.25 C------ end case "read topoFile"
193     ENDIF
194 jmc 1.24
195 jmc 1.25 C----- fill in the overlap :
196     _EXCH_XY_R4(Ro_surf, myThid )
197 jmc 1.24
198     c---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
199 adcroft 1.21
200 jmc 1.24 C------
201     C 3) Close the Domain (special configuration).
202     C------
203 adcroft 1.23 IF (groundAtK1) THEN
204 adcroft 1.21 DO bj = myByLo(myThid), myByHi(myThid)
205     DO bi = myBxLo(myThid), myBxHi(myThid)
206     DO j=1-Oly,sNy+Oly
207     DO i=1-Olx,sNx+Olx
208 jmc 1.24 iG = myXGlobalLo-1+(bi-1)*sNx+I
209     jG = myYGlobalLo-1+(bj-1)*sNy+J
210     C Test for eastern edge
211     c IF ( iG .EQ. Nx ) Ro_surf(i,j,bi,bj) = 0.
212     C Test for northern edge
213     c IF ( jG .EQ. Ny ) Ro_surf(i,j,bi,bj) = 0.
214 adcroft 1.23 IF (usingSphericalPolarGrid .AND. abs(yC(I,J,bi,bj)).GE.90. )
215 jmc 1.24 & Ro_surf(I,J,bi,bj) = rF(Nr+1)
216 adcroft 1.23 ENDDO
217     ENDDO
218     ENDDO
219     ENDDO
220     ELSE
221     DO bj = myByLo(myThid), myByHi(myThid)
222     DO bi = myBxLo(myThid), myBxHi(myThid)
223     DO j=1-Oly,sNy+Oly
224     DO i=1-Olx,sNx+Olx
225 jmc 1.24 iG = myXGlobalLo-1+(bi-1)*sNx+I
226     jG = myYGlobalLo-1+(bj-1)*sNy+J
227     C Test for eastern edge
228     c IF ( iG .EQ. Nx ) R_low(i,j,bi,bj) = 0.
229     C Test for northern edge
230     c IF ( jG .EQ. Ny ) R_low(i,j,bi,bj) = 0.
231 adcroft 1.23 IF (usingSphericalPolarGrid .AND. abs(yC(I,J,bi,bj)).GE.90. )
232 jmc 1.24 & R_low(I,J,bi,bj) = Ro_SeaLevel
233 adcroft 1.21 ENDDO
234     ENDDO
235     ENDDO
236     ENDDO
237     ENDIF
238 jmc 1.24
239     c _BEGIN_MASTER( myThid )
240     c CALL WRITE_FLD_XY_RS('Ro_surf',' ',Ro_surf,0,myThid)
241     c _END_MASTER(myThid)
242 adcroft 1.23
243 cnh 1.1 RETURN
244     END

  ViewVC Help
Powered by ViewVC 1.1.22