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

Annotation of /MITgcm/model/src/ini_spherical_polar_grid.F

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


Revision 1.2 - (hide annotations) (download)
Fri Apr 24 02:05:41 1998 UTC (26 years, 1 month ago) by cnh
Branch: MAIN
Changes since 1.1: +2 -2 lines
Further $Id to $Header conversions

1 cnh 1.2 C $Header: ini_spherical_polar_grid.F,v 1.1.1.1 1998/04/22 19:15:30 cnh Exp $
2 cnh 1.1
3     #include "CPP_EEOPTIONS.h"
4    
5     CStartOfInterface
6     SUBROUTINE INI_SPHERICAL_POLAR_GRID( myThid )
7     C /==========================================================\
8     C | SUBROUTINE INI_SPHERICAL_POLAR_GRID |
9     C | o Initialise model coordinate system |
10     C |==========================================================|
11     C | These arrays are used throughout the code in evaluating |
12     C | gradients, integrals and spatial avarages. This routine |
13     C | is called separately by each thread and initialise only |
14     C | the region of the domain it is "responsible" for. |
15     C | Notes: |
16     C | Two examples are included. One illustrates the |
17     C | initialisation of a cartesian grid. The other shows the |
18     C | inialisation of a spherical polar grid. Other orthonormal|
19     C | grids can be fitted into this design. In this case |
20     C | custom metric terms also need adding to account for the |
21     C | projections of velocity vectors onto these grids. |
22     C | The structure used here also makes it possible to |
23     C | implement less regular grid mappings. In particular |
24     C | o Schemes which leave out blocks of the domain that are |
25     C | all land could be supported. |
26     C | o Multi-level schemes such as icosohedral or cubic |
27     C | grid projections onto a sphere can also be fitted |
28     C | within the strategy we use. |
29     C | Both of the above also require modifying the support |
30     C | routines that map computational blocks to simulation |
31     C | domain blocks. |
32     C | Under the spherical polar grid mode primitive distances |
33     C | in X and Y are in degrees. Distance in Z are in m or Pa |
34     C | depending on the vertical gridding mode. |
35     C \==========================================================/
36    
37     C === Global variables ===
38     #include "SIZE.h"
39     #include "EEPARAMS.h"
40     #include "PARAMS.h"
41     #include "GRID.h"
42    
43     C == Routine arguments ==
44     C myThid - Number of this instance of INI_CARTESIAN_GRID
45     INTEGER myThid
46     CEndOfInterface
47    
48     C == Local variables ==
49     C xG, yG - Global coordinate location.
50     C zG
51     C xBase - South-west corner location for process.
52     C yBase
53     C zUpper - Work arrays for upper and lower
54     C zLower cell-face heights.
55     C phi - Temporary scalar
56     C iG, jG - Global coordinate index. Usually used to hold
57     C the south-west global coordinate of a tile.
58     C bi,bj - Loop counters
59     C zUpper - Temporary arrays holding z coordinates of
60     C zLower upper and lower faces.
61     C xBase - Lower coordinate for this threads cells
62     C yBase
63     C lat, latN, - Temporary variables used to hold latitude
64     C latS values.
65     C I,J,K
66     _RL xG, yG, zG
67     _RL phi
68     _RL zUpper(Nz), zLower(Nz)
69     _RL xBase, yBase
70     INTEGER iG, jG
71     INTEGER bi, bj
72     INTEGER I, J, K
73     _RL lat, latS, latN
74    
75     C-- Example of inialisation for spherical polar grid
76     C-- First set coordinates of cell centers
77     C This operation is only performed at start up so for more
78     C complex configurations it is usually OK to pass iG, jG to a custom
79     C function and have it return xG and yG.
80     C Set up my local grid first
81     C Note: In the spherical polar case delX and delY are given in
82     C degrees and are relative to some starting latitude and
83     C longitude - phiMin and thetaMin.
84     DO bj = myByLo(myThid), myByHi(myThid)
85     jG = myYGlobalLo + (bj-1)*sNy
86     DO bi = myBxLo(myThid), myBxHi(myThid)
87     iG = myXGlobalLo + (bi-1)*sNx
88     yBase = phiMin
89     xBase = thetaMin
90     DO i=1,iG-1
91     xBase = xBase + delX(i)
92     ENDDO
93     DO j=1,jG-1
94     yBase = yBase + delY(j)
95     ENDDO
96     yG = yBase
97     DO J=1,sNy
98     xG = xBase
99     DO I=1,sNx
100     xc(I,J,bi,bj) = xG + delX(iG+i-1)*0.5 _d 0
101     yc(I,J,bi,bj) = yG + delY(jG+j-1)*0.5 _d 0
102     xG = xG + delX(iG+I-1)
103     dxF(I,J,bi,bj) = delX(iG+i-1)*deg2rad*rSphere*COS(yc(I,J,bi,bj)*deg2rad)
104     dyF(I,J,bi,bj) = delY(jG+j-1)*deg2rad*rSphere
105     ENDDO
106     yG = yG + delY(jG+J-1)
107     ENDDO
108     ENDDO
109     ENDDO
110     C Now sync. and get edge regions from other threads and/or processes.
111     C Note: We could just set the overlap regions ourselves here but
112     C exchanging edges is safer and is good practice!
113     _EXCH_XY_R4( xc, myThid )
114     _EXCH_XY_R4( yc, myThid )
115     _EXCH_XY_R4(dxF, myThid )
116     _EXCH_XY_R4(dyF, myThid )
117    
118     C-- Calculate separation between other points
119     C dxG, dyG are separations between cell corners along cell faces.
120     DO bj = myByLo(myThid), myByHi(myThid)
121     DO bi = myBxLo(myThid), myBxHi(myThid)
122     DO J=1,sNy
123     DO I=1,sNx
124     jG = myYGlobalLo + (bj-1)*sNy + J-1
125     iG = myXGlobalLo + (bi-1)*sNx + I-1
126     lat = yc(I,J,bi,bj)-delY(jG) * 0.5 _d 0
127     dxG(I,J,bi,bj) = rSphere*COS(lat*deg2rad)*delX(iG)*deg2rad
128     dyG(I,J,bi,bj) = (dyF(I,J,bi,bj)+dyF(I-1,J,bi,bj))*0.5 _d 0
129     ENDDO
130     ENDDO
131     ENDDO
132     ENDDO
133     _EXCH_XY_R4(dxG, myThid )
134     _EXCH_XY_R4(dyG, myThid )
135     C dxV, dyU are separations between velocity points along cell faces.
136     DO bj = myByLo(myThid), myByHi(myThid)
137     DO bi = myBxLo(myThid), myBxHi(myThid)
138     DO J=1,sNy
139     DO I=1,sNx
140     dxV(I,J,bi,bj) = (dxG(I,J,bi,bj)+dxG(I-1,J,bi,bj))*0.5 _d 0
141     dyU(I,J,bi,bj) = (dyG(I,J,bi,bj)+dyG(I,J-1,bi,bj))*0.5 _d 0
142     ENDDO
143     ENDDO
144     ENDDO
145     ENDDO
146     _EXCH_XY_R4(dxV, myThid )
147     _EXCH_XY_R4(dyU, myThid )
148     C dxC, dyC is separation between cell centers
149     DO bj = myByLo(myThid), myByHi(myThid)
150     DO bi = myBxLo(myThid), myBxHi(myThid)
151     DO J=1,sNy
152     DO I=1,sNx
153     dxC(I,J,bi,bj) = (dxF(I,J,bi,bj)+dxF(I-1,J,bi,bj))*0.5 _d 0
154     dyC(I,J,bi,bj) = (dyF(I,J,bi,bj)+dyF(I,J-1,bi,bj))*0.5 _d 0
155     ENDDO
156     ENDDO
157     ENDDO
158     ENDDO
159     _EXCH_XY_R4(dxC, myThid )
160     _EXCH_XY_R4(dyC, myThid )
161     C Calculate recipricols
162     DO bj = myByLo(myThid), myByHi(myThid)
163     DO bi = myBxLo(myThid), myBxHi(myThid)
164     DO J=1,sNy
165     DO I=1,sNx
166     rDxG(I,J,bi,bj)=1.d0/dxG(I,J,bi,bj)
167     rDyG(I,J,bi,bj)=1.d0/dyG(I,J,bi,bj)
168     rDxC(I,J,bi,bj)=1.d0/dxC(I,J,bi,bj)
169     rDyC(I,J,bi,bj)=1.d0/dyC(I,J,bi,bj)
170     rDxF(I,J,bi,bj)=1.d0/dxF(I,J,bi,bj)
171     rDyF(I,J,bi,bj)=1.d0/dyF(I,J,bi,bj)
172     rDxV(I,J,bi,bj)=1.d0/dxV(I,J,bi,bj)
173     rDyU(I,J,bi,bj)=1.d0/dyU(I,J,bi,bj)
174     ENDDO
175     ENDDO
176     ENDDO
177     ENDDO
178     _EXCH_XY_R4(rDxG, myThid )
179     _EXCH_XY_R4(rDyG, myThid )
180     _EXCH_XY_R4(rDxC, myThid )
181     _EXCH_XY_R4(rDyC, myThid )
182     _EXCH_XY_R4(rDxF, myThid )
183     _EXCH_XY_R4(rDyF, myThid )
184     _EXCH_XY_R4(rDxV, myThid )
185     _EXCH_XY_R4(rDyU, myThid )
186     C Calculate vertical face area
187     DO bj = myByLo(myThid), myByHi(myThid)
188     DO bi = myBxLo(myThid), myBxHi(myThid)
189     DO J=1,sNy
190     DO I=1,sNx
191     jG = myYGlobalLo + (bj-1)*sNy + J-1
192     latS = yc(i,j,bi,bj)-delY(jG)*0.5 _d 0
193     latN = yc(i,j,bi,bj)+delY(jG)*0.5 _d 0
194     zA(I,J,bi,bj) = dyF(I,J,bi,bj)
195     & *rSphere*(SIN(latN*deg2rad)-SIN(latS*deg2rad))
196     ENDDO
197     ENDDO
198     ENDDO
199     ENDDO
200     C
201     RETURN
202     END
203    
204 cnh 1.2 C $Id: ini_spherical_polar_grid.F,v 1.1.1.1 1998/04/22 19:15:30 cnh Exp $

  ViewVC Help
Powered by ViewVC 1.1.22