/[MITgcm]/MITgcm/verification/aim.5l_Equatorial_Channel/code/ini_depths.F
ViewVC logotype

Diff of /MITgcm/verification/aim.5l_Equatorial_Channel/code/ini_depths.F

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

revision 1.3 by cnh, Sun Feb 4 14:38:51 2001 UTC revision 1.4 by adcroft, Tue May 29 14:01:46 2001 UTC
# Line 27  C     === Global variables === Line 27  C     === Global variables ===
27  #include "EEPARAMS.h"  #include "EEPARAMS.h"
28  #include "PARAMS.h"  #include "PARAMS.h"
29  #include "GRID.h"  #include "GRID.h"
30    #include "INI_DEPTHS.h"
31    
32  C     == Routine arguments ==  C     == Routine arguments ==
33  C     myThid -  Number of this instance of INI_DEPTHS  C     myThid -  Number of this instance of INI_DEPTHS
# Line 37  C     == Local variables == Line 38  C     == Local variables ==
38  C     iG, jG - Global coordinate index  C     iG, jG - Global coordinate index
39  C     bi,bj  - Loop counters  C     bi,bj  - Loop counters
40  C     I,J,K  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)  C     Hdefault - default r-coordinate of the lower boundary (=ground)
43  C                (=minus(Total_depth) in the ocean model)  C                (=minus(Total_depth) in the ocean model)
44  C                (=Total Pressure at Sea Level in the atmos model)  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  C     oldPrec - Temporary used in controlling binary input dataset precision
50        INTEGER iG, jG        INTEGER iG, jG
51        INTEGER bi, bj        INTEGER bi, bj
52        INTEGER  I, J, K        INTEGER  I, J
53        _RL Hdefault        _RL Hdefault
54    
55        _BARRIER        _BARRIER
# Line 52  C      Set up a flat bottom box with dou Line 58  C      Set up a flat bottom box with dou
58  C      H is the basic variable from which other terms are derived. It  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  C      is the term that would be set from an external file for a
60  C      realistic problem.  C      realistic problem.
61         Hdefault = rF(1)         IF (groundAtK1) THEN
62         IF (.NOT. groundAtK1) THEN           Hdefault = Ro_SeaLevel
63          DO K = 1, Nr         ELSE
64           Hdefault = Hdefault - rkfac*drF(K)           Hdefault = rF(Nr+1)
         ENDDO  
65         ENDIF         ENDIF
66         DO bj = myByLo(myThid), myByHi(myThid)         DO bj = myByLo(myThid), myByHi(myThid)
67          DO bi = myBxLo(myThid), myBxHi(myThid)          DO bi = myBxLo(myThid), myBxHi(myThid)
# Line 69  C          Default depth of full domain Line 74  C          Default depth of full domain
74  C          Test for eastern edge  C          Test for eastern edge
75  C          IF ( iG .EQ. nX ) H(i,j,bi,bj) = 0.  C          IF ( iG .EQ. nX ) H(i,j,bi,bj) = 0.
76  C          Test for northern edge  C          Test for northern edge
77             IF ( jG .EQ. nY ) H(i,j,bi,bj) = 0.  C          IF ( jG .EQ. nY ) H(i,j,bi,bj) = 0.
 C          H(i,j,bi,bj) = 0.  
78            ENDDO            ENDDO
79           ENDDO           ENDDO
80          ENDDO          ENDDO
# Line 91  C    &                    'RS', 1, H, 1, Line 95  C    &                    'RS', 1, H, 1,
95    
96        _EXCH_XY_R4(    H, myThid )        _EXCH_XY_R4(    H, myThid )
97    
98        IF (usingSphericalPolarGrid) THEN  C
99    C     CALL PLOT_FIELD_XYRS(H,'Model depths (ini_depths)',1,myThid)
100    C
101          IF (groundAtK1) THEN
102         DO bj = myByLo(myThid), myByHi(myThid)         DO bj = myByLo(myThid), myByHi(myThid)
103          DO bi = myBxLo(myThid), myBxHi(myThid)          DO bi = myBxLo(myThid), myBxHi(myThid)
104           DO j=1-Oly,sNy+Oly           DO j=1-Oly,sNy+Oly
105            DO i=1-Olx,sNx+Olx            DO i=1-Olx,sNx+Olx
106             IF (abs(yC(I,J,bi,bj)).GE.90.) H(I,J,bi,bj)=0.             R_low(I,J,bi,bj) = rF(Nr+1)
107               Ro_surf(I,J,bi,bj) = H(I,J,bi,bj)
108    c- Domain : Symetric % Eq. & closed at N & S boundaries:
109               IF (usingSphericalPolarGrid .AND.
110         &         abs(yC(I,J,bi,bj)).GE.-phiMin) Ro_surf(I,J,bi,bj) = 0.
111               IF (usingSphericalPolarGrid .AND. abs(yC(I,J,bi,bj)).GE.90. )
112         &     THEN
113                Ro_surf(I,J,bi,bj) = 0.
114               ENDIF
115              ENDDO
116             ENDDO
117            ENDDO
118           ENDDO
119          ELSE
120           DO bj = myByLo(myThid), myByHi(myThid)
121            DO bi = myBxLo(myThid), myBxHi(myThid)
122             DO j=1-Oly,sNy+Oly
123              DO i=1-Olx,sNx+Olx
124               R_low(I,J,bi,bj) = H(I,J,bi,bj)
125               Ro_surf(I,J,bi,bj) = Ro_SeaLevel
126    c- Domain : Symetric % Eq. & closed at N & S boundaries:
127               IF (usingSphericalPolarGrid .AND.
128         &         abs(yC(I,J,bi,bj)).GE.-phiMin) R_low(I,J,bi,bj) = 0.
129               IF (usingSphericalPolarGrid .AND. abs(yC(I,J,bi,bj)).GE.90. )
130         &     THEN
131                R_low(I,J,bi,bj) = 0.
132               ENDIF
133            ENDDO            ENDDO
134           ENDDO           ENDDO
135          ENDDO          ENDDO
136         ENDDO         ENDDO
137        ENDIF        ENDIF
138  C  
       CALL PLOT_FIELD_XYRS(H,'Model depths (ini_depths)',1,myThid)  
 C  
139        RETURN        RETURN
140        END        END

Legend:
Removed from v.1.3  
changed lines
  Added in v.1.4

  ViewVC Help
Powered by ViewVC 1.1.22