/[MITgcm]/MITgcm/model/src/ini_depths.F
ViewVC logotype

Contents of /MITgcm/model/src/ini_depths.F

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


Revision 1.10 - (show annotations) (download)
Wed Jul 15 22:20:02 1998 UTC (25 years, 10 months ago) by adcroft
Branch: MAIN
CVS Tags: branch-point-rdot
Branch point for: branch-rdot
Changes since 1.9: +11 -1 lines
This is a bit confusing. We use Z coordinates with Z=0 at the top
and z=H at the bottom.  The depth of the model (internally) is positive
but we wanted to use negative depths for input files. That way we can
use the same bathymetry file for the atmospheric model which would have
positive elevation for the bottom boundary. Oh boy!
I've added a sign change immediately after the read of bathyFile
to make this work but we ought to write this all down somewhere...

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_depths.F,v 1.9 1998/07/02 15:46:21 adcroft Exp $
2
3 #include "CPP_EEOPTIONS.h"
4
5 CStartOfInterface
6 SUBROUTINE INI_DEPTHS( myThid )
7 C /==========================================================\
8 C | SUBROUTINE INI_DEPTHS |
9 C | o Initialise map of model depths |
10 C |==========================================================|
11 C | The depths of the bottom of the model is specified in |
12 C | terms of an XY map with one depth for each column of |
13 C | grid cells. Depths do not have to coincide with the |
14 C | model levels. The model's lopping algorithm makes it |
15 C | possible to represent arbitrary depths. |
16 C | The mode depths map also influences the models topology |
17 C | By default the model domain wraps around in X and Y. |
18 C | This default doubly periodic topology is "supressed" |
19 C | if a depth map is defined which closes off all wrap |
20 C | around flow. |
21 C \==========================================================/
22
23 C === Global variables ===
24 #include "SIZE.h"
25 #include "EEPARAMS.h"
26 #include "PARAMS.h"
27 #include "GRID.h"
28
29 C == Routine arguments ==
30 C myThid - Number of this instance of INI_DEPTHS
31 INTEGER myThid
32 CEndOfInterface
33
34 C == Local variables ==
35 C iG, jG - Global coordinate index
36 C bi,bj - Loop counters
37 C I,J,K
38 C phi - total depth of model
39 INTEGER iG, jG
40 INTEGER bi, bj
41 INTEGER I, J, K
42 _RL phi
43
44 _BARRIER
45 IF ( bathyFile .EQ. ' ' ) THEN
46 C Set up a flat bottom box with doubly periodic topology.
47 C H is the basic variable from which other terms are derived. It
48 C is the term that would be set from an external file for a
49 C realistic problem.
50 phi = zFace(nZ+1)
51 DO bj = myByLo(myThid), myByHi(myThid)
52 DO bi = myBxLo(myThid), myBxHi(myThid)
53 DO j=1,sNy
54 DO i=1,sNx
55 iG = myXGlobalLo-1+(bi-1)*sNx+I
56 jG = myYGlobalLo-1+(bj-1)*sNy+J
57 C Default depth of full domain
58 H(i,j,bi,bj) = phi
59 C Test for eastern edge
60 IF ( iG .EQ. nX ) H(i,j,bi,bj) = 0.
61 C Test for northern edge
62 IF ( jG .EQ. nY ) H(i,j,bi,bj) = 0.
63 C Island
64 IF ( iG .EQ. 1 .AND.
65 & jG .EQ. 24 ) H(i,j,bi,bj) = 0.75*phi
66 ENDDO
67 ENDDO
68 ENDDO
69 ENDDO
70 ELSE
71 _BEGIN_MASTER( myThid )
72 CALL READ_FLD_XY_RS( bathyFile, ' ', H, 0, myThid )
73 _END_MASTER(myThid)
74 C-- Reverse sign if necessary (We assume H>0)
75 DO bj = myByLo(myThid), myByHi(myThid)
76 DO bi = myBxLo(myThid), myBxHi(myThid)
77 DO j=1,sNy
78 DO i=1,sNx
79 H(i,j,bi,bj)=-H(i,j,bi,bj)
80 ENDDO
81 ENDDO
82 ENDDO
83 ENDDO
84 ENDIF
85 _EXCH_XY_R4( H, myThid )
86 C
87 CALL PLOT_FIELD_XYRS( H, 'Model depths' , 1, myThid )
88 C
89 RETURN
90 END

  ViewVC Help
Powered by ViewVC 1.1.22