/[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.1 - (hide annotations) (download)
Fri Sep 28 02:28:10 2001 UTC (22 years, 7 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint46b_post, checkpoint46k_post, checkpoint44h_pre, release1_p12, release1_p13, release1_p10, release1_p16, release1_p15, release1_p11, checkpoint46j_post, checkpoint44e_post, checkpoint44f_pre, checkpoint47a_post, checkpoint46c_post, checkpoint46f_post, checkpoint46l_pre, checkpoint46a_post, checkpoint46n_post, release1_beta1, checkpoint46d_pre, checkpoint46e_post, release1-branch_tutorials, release1_p14, checkpoint44g_post, checkpoint46h_pre, checkpoint44h_post, checkpoint46l_post, checkpoint46e_pre, checkpoint43a-release1mods, checkpoint45d_post, checkpoint46j_pre, checkpoint45b_post, checkpoint46b_pre, chkpt44a_pre, release1-branch-end, release1_final_v1, checkpoint46, checkpoint44, checkpoint45, checkpoint44f_post, checkpoint47b_post, release1_p17, release1_b1, checkpoint44b_post, chkpt44d_post, release1_p8, release1_p9, checkpoint46m_post, checkpoint46g_pre, release1_p2, release1_p3, release1_p4, release1_p6, chkpt44a_post, checkpoint44b_pre, release1_p1, checkpoint46a_pre, checkpoint45c_post, release1_p5, checkpoint44e_pre, chkpt44c_pre, release1_p7, checkpoint46d_post, checkpoint46g_post, release1_p13_pre, release1_p12_pre, checkpoint45a_post, checkpoint46c_pre, checkpoint43, release1-branch_branchpoint, checkpoint46i_post, checkpoint47, chkpt44c_post, checkpoint46h_post, release1_chkpt44d_post
Branch point for: release1_coupled, release1_final, release1-branch, release1, release1_50yr
Adding test of advection schemes with lopped cells. 2-D x-z

1 adcroft 1.1 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    
4     #include "CPP_OPTIONS.h"
5    
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     C | model levels. The model lopping algorithm makes it |
16     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     IMPLICIT NONE
24    
25     C === Global variables ===
26     #include "SIZE.h"
27     #include "EEPARAMS.h"
28     #include "PARAMS.h"
29     #include "GRID.h"
30     #include "INI_DEPTHS.h"
31    
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     C I,J,K
41     C H - Depth of base of fluid from upper surface f[X,Y] (m).
42     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     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     C oldPrec - Temporary used in controlling binary input dataset precision
50     INTEGER iG, jG
51     INTEGER bi, bj
52     INTEGER I, J
53     _RL Hdefault
54    
55     _BARRIER
56     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     IF (groundAtK1) THEN
62     Hdefault = Ro_SeaLevel
63     ELSE
64     Hdefault = rF(Nr+1)
65     ENDIF
66     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     H(i,j,bi,bj) = Hdefault
74     & *(1.-0.5*XC(I,J,bi,bj)/(float(Nx)*DelX(1)))
75     C Test for eastern edge
76     C IF ( iG .EQ. nX ) H(i,j,bi,bj) = 0.
77     C Test for northern edge
78     C IF ( jG .EQ. nY ) H(i,j,bi,bj) = 0.
79     ENDDO
80     ENDDO
81     ENDDO
82     ENDDO
83     ELSE
84     _BEGIN_MASTER( myThid )
85     C Read the bathymetry using the mid-level I/O pacakage read_write_rec
86     C The 0 is the "iteration" argument. The 1 is the record number.
87     CALL READ_REC_XY_RS( bathyFile, H, 1, 0, myThid )
88     C Read the bathymetry using the mid-level I/O pacakage read_write_fld
89     C The 0 is the "iteration" argument. The ' ' is an empty suffix
90     C CALL READ_FLD_XY_RS( bathyFile, ' ', H, 0, myThid )
91     C Read the bathymetry using the low-level I/O package
92     C CALL MDSREADFIELD( bathyFile, readBinaryPrec,
93     C & 'RS', 1, H, 1, myThid )
94     _END_MASTER(myThid)
95     ENDIF
96    
97     _EXCH_XY_R4( H, myThid )
98    
99     C
100     C CALL PLOT_FIELD_XYRS(H,'Model depths (ini_depths)',1,myThid)
101     C
102     IF (groundAtK1) THEN
103     DO bj = myByLo(myThid), myByHi(myThid)
104     DO bi = myBxLo(myThid), myBxHi(myThid)
105     DO j=1-Oly,sNy+Oly
106     DO i=1-Olx,sNx+Olx
107     R_low(I,J,bi,bj) = rF(Nr+1)
108     Ro_surf(I,J,bi,bj) = H(I,J,bi,bj)
109     IF (usingSphericalPolarGrid .AND. abs(yC(I,J,bi,bj)).GE.90. )
110     & THEN
111     Ro_surf(I,J,bi,bj) = 0.
112     ENDIF
113     ENDDO
114     ENDDO
115     ENDDO
116     ENDDO
117     ELSE
118     DO bj = myByLo(myThid), myByHi(myThid)
119     DO bi = myBxLo(myThid), myBxHi(myThid)
120     DO j=1-Oly,sNy+Oly
121     DO i=1-Olx,sNx+Olx
122     R_low(I,J,bi,bj) = H(I,J,bi,bj)
123     Ro_surf(I,J,bi,bj) = Ro_SeaLevel
124     IF (usingSphericalPolarGrid .AND. abs(yC(I,J,bi,bj)).GE.90. )
125     & THEN
126     R_low(I,J,bi,bj) = 0.
127     ENDIF
128     ENDDO
129     ENDDO
130     ENDDO
131     ENDDO
132     ENDIF
133    
134     RETURN
135     END

  ViewVC Help
Powered by ViewVC 1.1.22