/[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.22 - (hide annotations) (download)
Sun Feb 4 14:38:47 2001 UTC (23 years, 4 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint38, c37_adj, checkpoint39, checkpoint37, checkpoint36, checkpoint35
Branch point for: pre38
Changes since 1.21: +2 -1 lines
Made sure each .F and .h file had
the CVS keywords Header and Name at its start.
Most had header but very few currently have Name, so
lots of changes!

1 cnh 1.22 C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_depths.F,v 1.21 2001/02/02 21:04:48 adcroft 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    
31     C == Routine arguments ==
32     C myThid - Number of this instance of INI_DEPTHS
33     INTEGER myThid
34     CEndOfInterface
35    
36     C == Local variables ==
37     C iG, jG - Global coordinate index
38     C bi,bj - Loop counters
39 adcroft 1.21 C I,J,K
40     C Hdefault - default r-coordinate of the lower boundary (=ground)
41     C (=minus(Total_depth) in the ocean model)
42     C (=Total Pressure at Sea Level in the atmos model)
43 cnh 1.13 C oldPrec - Temporary used in controlling binary input dataset precision
44 cnh 1.1 INTEGER iG, jG
45     INTEGER bi, bj
46 adcroft 1.21 INTEGER I, J, K
47     _RL Hdefault
48 cnh 1.1
49     _BARRIER
50 cnh 1.4 IF ( bathyFile .EQ. ' ' ) THEN
51     C Set up a flat bottom box with doubly periodic topology.
52     C H is the basic variable from which other terms are derived. It
53     C is the term that would be set from an external file for a
54     C realistic problem.
55 adcroft 1.21 Hdefault = rF(1)
56     IF (.NOT. groundAtK1) THEN
57     DO K = 1, Nr
58     Hdefault = Hdefault - rkfac*drF(K)
59     ENDDO
60     ENDIF
61 cnh 1.4 DO bj = myByLo(myThid), myByHi(myThid)
62     DO bi = myBxLo(myThid), myBxHi(myThid)
63     DO j=1,sNy
64     DO i=1,sNx
65     iG = myXGlobalLo-1+(bi-1)*sNx+I
66     jG = myYGlobalLo-1+(bj-1)*sNy+J
67     C Default depth of full domain
68 adcroft 1.21 H(i,j,bi,bj) = Hdefault
69 cnh 1.4 C Test for eastern edge
70 adcroft 1.18 C IF ( iG .EQ. nX ) H(i,j,bi,bj) = 0.
71 cnh 1.4 C Test for northern edge
72 adcroft 1.18 C IF ( jG .EQ. nY ) H(i,j,bi,bj) = 0.
73 cnh 1.4 ENDDO
74 cnh 1.1 ENDDO
75     ENDDO
76     ENDDO
77 cnh 1.4 ELSE
78     _BEGIN_MASTER( myThid )
79 adcroft 1.19 C Read the bathymetry using the mid-level I/O pacakage read_write_rec
80     C The 0 is the "iteration" argument. The 1 is the record number.
81     CALL READ_REC_XY_RS( bathyFile, H, 1, 0, myThid )
82     C Read the bathymetry using the mid-level I/O pacakage read_write_fld
83     C The 0 is the "iteration" argument. The ' ' is an empty suffix
84     C CALL READ_FLD_XY_RS( bathyFile, ' ', H, 0, myThid )
85     C Read the bathymetry using the low-level I/O package
86     C CALL MDSREADFIELD( bathyFile, readBinaryPrec,
87     C & 'RS', 1, H, 1, myThid )
88 cnh 1.4 _END_MASTER(myThid)
89     ENDIF
90 adcroft 1.21
91 cnh 1.4 _EXCH_XY_R4( H, myThid )
92 adcroft 1.21
93     IF (usingSphericalPolarGrid) THEN
94     DO bj = myByLo(myThid), myByHi(myThid)
95     DO bi = myBxLo(myThid), myBxHi(myThid)
96     DO j=1-Oly,sNy+Oly
97     DO i=1-Olx,sNx+Olx
98     IF (abs(yC(I,J,bi,bj)).GE.90.) H(I,J,bi,bj)=0.
99     ENDDO
100     ENDDO
101     ENDDO
102     ENDDO
103     ENDIF
104 cnh 1.5 C
105 adcroft 1.21 CALL PLOT_FIELD_XYRS(H,'Model depths (ini_depths)',1,myThid)
106 cnh 1.1 C
107     RETURN
108     END

  ViewVC Help
Powered by ViewVC 1.1.22