/[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.23 - (hide annotations) (download)
Tue May 29 14:01:37 2001 UTC (23 years ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint40pre1
Changes since 1.22: +39 -13 lines
Merge from branch pre38:
 o essential mods for cubed sphere
 o debugged atmosphere, dynamcis + physics (aim)
 o new packages (mom_vecinv, mom_fluxform, ...)

1 adcroft 1.23 C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_depths.F,v 1.22.2.5 2001/04/13 10:55:23 cnh Exp $
2     C $Name: $
3 cnh 1.1
4 cnh 1.16 #include "CPP_OPTIONS.h"
5 cnh 1.1
6     CStartOfInterface
7     SUBROUTINE INI_DEPTHS( myThid )
8     C /==========================================================\
9     C | SUBROUTINE INI_DEPTHS |
10     C | o Initialise map of model depths |
11     C |==========================================================|
12     C | The depths of the bottom of the model is specified in |
13     C | terms of an XY map with one depth for each column of |
14     C | grid cells. Depths do not have to coincide with the |
15 cnh 1.15 C | model levels. The model lopping algorithm makes it |
16 cnh 1.1 C | possible to represent arbitrary depths. |
17     C | The mode depths map also influences the models topology |
18     C | By default the model domain wraps around in X and Y. |
19     C | This default doubly periodic topology is "supressed" |
20     C | if a depth map is defined which closes off all wrap |
21     C | around flow. |
22     C \==========================================================/
23 adcroft 1.17 IMPLICIT NONE
24 cnh 1.1
25     C === Global variables ===
26     #include "SIZE.h"
27     #include "EEPARAMS.h"
28     #include "PARAMS.h"
29     #include "GRID.h"
30 adcroft 1.23 #include "INI_DEPTHS.h"
31 cnh 1.1
32     C == Routine arguments ==
33     C myThid - Number of this instance of INI_DEPTHS
34     INTEGER myThid
35     CEndOfInterface
36    
37     C == Local variables ==
38     C iG, jG - Global coordinate index
39     C bi,bj - Loop counters
40 adcroft 1.21 C I,J,K
41 adcroft 1.23 C H - Depth of base of fluid from upper surface f[X,Y] (m).
42 adcroft 1.21 C Hdefault - default r-coordinate of the lower boundary (=ground)
43     C (=minus(Total_depth) in the ocean model)
44     C (=Total Pressure at Sea Level in the atmos model)
45 adcroft 1.23 C--------------------
46     C NOTE: will change soon: 2 separed files for r_lower and r_surface boudaries
47     C and for the atmosphere, topography will be defined in term of height
48     C--------------------
49 cnh 1.13 C oldPrec - Temporary used in controlling binary input dataset precision
50 cnh 1.1 INTEGER iG, jG
51     INTEGER bi, bj
52 adcroft 1.23 INTEGER I, J
53 adcroft 1.21 _RL Hdefault
54 cnh 1.1
55     _BARRIER
56 cnh 1.4 IF ( bathyFile .EQ. ' ' ) THEN
57     C Set up a flat bottom box with doubly periodic topology.
58     C H is the basic variable from which other terms are derived. It
59     C is the term that would be set from an external file for a
60     C realistic problem.
61 adcroft 1.23 IF (groundAtK1) THEN
62     Hdefault = Ro_SeaLevel
63     ELSE
64     Hdefault = rF(Nr+1)
65 adcroft 1.21 ENDIF
66 cnh 1.4 DO bj = myByLo(myThid), myByHi(myThid)
67     DO bi = myBxLo(myThid), myBxHi(myThid)
68     DO j=1,sNy
69     DO i=1,sNx
70     iG = myXGlobalLo-1+(bi-1)*sNx+I
71     jG = myYGlobalLo-1+(bj-1)*sNy+J
72     C Default depth of full domain
73 adcroft 1.21 H(i,j,bi,bj) = Hdefault
74 cnh 1.4 C Test for eastern edge
75 adcroft 1.18 C IF ( iG .EQ. nX ) H(i,j,bi,bj) = 0.
76 cnh 1.4 C Test for northern edge
77 adcroft 1.18 C IF ( jG .EQ. nY ) H(i,j,bi,bj) = 0.
78 cnh 1.4 ENDDO
79 cnh 1.1 ENDDO
80     ENDDO
81     ENDDO
82 cnh 1.4 ELSE
83     _BEGIN_MASTER( myThid )
84 adcroft 1.19 C Read the bathymetry using the mid-level I/O pacakage read_write_rec
85     C The 0 is the "iteration" argument. The 1 is the record number.
86     CALL READ_REC_XY_RS( bathyFile, H, 1, 0, myThid )
87     C Read the bathymetry using the mid-level I/O pacakage read_write_fld
88     C The 0 is the "iteration" argument. The ' ' is an empty suffix
89     C CALL READ_FLD_XY_RS( bathyFile, ' ', H, 0, myThid )
90     C Read the bathymetry using the low-level I/O package
91     C CALL MDSREADFIELD( bathyFile, readBinaryPrec,
92     C & 'RS', 1, H, 1, myThid )
93 cnh 1.4 _END_MASTER(myThid)
94     ENDIF
95 adcroft 1.21
96 cnh 1.4 _EXCH_XY_R4( H, myThid )
97 adcroft 1.21
98 adcroft 1.23 C
99     C CALL PLOT_FIELD_XYRS(H,'Model depths (ini_depths)',1,myThid)
100     C
101     IF (groundAtK1) THEN
102 adcroft 1.21 DO bj = myByLo(myThid), myByHi(myThid)
103     DO bi = myBxLo(myThid), myBxHi(myThid)
104     DO j=1-Oly,sNy+Oly
105     DO i=1-Olx,sNx+Olx
106 adcroft 1.23 R_low(I,J,bi,bj) = rF(Nr+1)
107     Ro_surf(I,J,bi,bj) = H(I,J,bi,bj)
108     IF (usingSphericalPolarGrid .AND. abs(yC(I,J,bi,bj)).GE.90. )
109     & THEN
110     Ro_surf(I,J,bi,bj) = 0.
111     ENDIF
112     ENDDO
113     ENDDO
114     ENDDO
115     ENDDO
116     ELSE
117     DO bj = myByLo(myThid), myByHi(myThid)
118     DO bi = myBxLo(myThid), myBxHi(myThid)
119     DO j=1-Oly,sNy+Oly
120     DO i=1-Olx,sNx+Olx
121     R_low(I,J,bi,bj) = H(I,J,bi,bj)
122     Ro_surf(I,J,bi,bj) = Ro_SeaLevel
123     IF (usingSphericalPolarGrid .AND. abs(yC(I,J,bi,bj)).GE.90. )
124     & THEN
125     R_low(I,J,bi,bj) = 0.
126     ENDIF
127 adcroft 1.21 ENDDO
128     ENDDO
129     ENDDO
130     ENDDO
131     ENDIF
132 adcroft 1.23
133 cnh 1.1 RETURN
134     END

  ViewVC Help
Powered by ViewVC 1.1.22