1 |
C $Header: /u/u3/gcmpack/MITgcm/model/src/ini_pressure.F,v 1.5.2.1 2003/10/02 18:10:45 edhill Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: INI_PRESSURE |
8 |
C !INTERFACE: |
9 |
SUBROUTINE INI_PRESSURE( myThid ) |
10 |
C !DESCRIPTION: \bv |
11 |
C *==========================================================* |
12 |
C | SUBROUTINE INI_PRESSURE |
13 |
C | o initialise the pressure field consistently with |
14 |
C | temperature and salinity |
15 |
C | - needs to be called after ini_theta, ini_salt, and |
16 |
C | ini_psurf |
17 |
C | - does not include surface pressure loading, because |
18 |
C | that is only available until after |
19 |
C | CALL packages_init_variables |
20 |
C *==========================================================* |
21 |
C \ev |
22 |
|
23 |
C !USES: |
24 |
IMPLICIT NONE |
25 |
C == Global variables == |
26 |
#include "SIZE.h" |
27 |
#include "EEPARAMS.h" |
28 |
#include "PARAMS.h" |
29 |
#include "EOS.h" |
30 |
#include "GRID.h" |
31 |
#include "DYNVARS.h" |
32 |
|
33 |
C !INPUT/OUTPUT PARAMETERS: |
34 |
C == Routine arguments == |
35 |
C myThid - Number of this instance of INI_PRESSURE |
36 |
INTEGER myThid |
37 |
|
38 |
C !LOCAL VARIABLES: |
39 |
C == Local variables == |
40 |
C dPhiHydX,Y :: Gradient (X & Y directions) of Hyd. Potential |
41 |
C bi,bj - Loop counters |
42 |
C I,J,K |
43 |
INTEGER bi, bj |
44 |
INTEGER I, J, K |
45 |
INTEGER iMin, iMax, jMin, jMax, npiter |
46 |
_RL PhiHydF (1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
47 |
_RL PhiHydC (1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
48 |
_RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
49 |
_RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
50 |
_RL oldPhi (1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
51 |
_RL count, rmspp, rmsppold |
52 |
|
53 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
54 |
CEOP |
55 |
iMin = 1-OLx |
56 |
iMax = sNx+OLx |
57 |
jMin = 1-OLy |
58 |
jMax = sNy+OLy |
59 |
|
60 |
WRITE(msgBuf,'(a)') |
61 |
& 'Start initial hydrostatic pressure computation' |
62 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
63 |
& SQUEEZE_RIGHT , 1) |
64 |
DO bj = myByLo(myThid), myByHi(myThid) |
65 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
66 |
DO k = 1,Nr |
67 |
DO j=1-OLy,sNy+OLy |
68 |
DO i=1-OLx,sNx+OLx |
69 |
totPhiHyd(i,j,k,bi,bj) = 0. _d 0 |
70 |
ENDDO |
71 |
ENDDO |
72 |
ENDDO |
73 |
ENDDO |
74 |
ENDDO |
75 |
|
76 |
c IF ( startTime .NE. 0. ) RETURN |
77 |
c IF ( buoyancyRelation .EQ. 'OCEANIC' ) THEN |
78 |
IF ( useDynP_inEos_Zc ) THEN |
79 |
|
80 |
#ifndef ALLOW_AUTODIFF_TAMC |
81 |
cph-- deal with this iterative loop for AD once it will |
82 |
cph-- really be needed; |
83 |
cph-- would need storing of totPhiHyd for each npiter |
84 |
|
85 |
rmspp = 1. _d 0 |
86 |
rmsppold = 0. _d 0 |
87 |
npiter = 0 |
88 |
|
89 |
C iterate pressure p = integral of (g*rho(p)*dz) |
90 |
DO npiter= 1, 15 |
91 |
rmsppold = rmspp |
92 |
rmspp = 0. _d 0 |
93 |
count = 0. |
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 |
phiHydF(i,j) = 0. _d 0 |
99 |
ENDDO |
100 |
ENDDO |
101 |
DO k = 1, Nr |
102 |
C for each level save old pressure and compute new pressure |
103 |
DO j=jMin,jMax |
104 |
DO i=iMin,iMax |
105 |
oldPhi(i,j) = totPhiHyd(i,j,k,bi,bj) |
106 |
ENDDO |
107 |
ENDDO |
108 |
CALL CALC_PHI_HYD( |
109 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
110 |
I theta, salt, |
111 |
U phiHydF, |
112 |
O phiHydC, dPhiHydX, dPhiHydY, |
113 |
I startTime, nIter0, myThid) |
114 |
C compute convergence criterion |
115 |
DO j=jMin,jMax |
116 |
DO i=iMin,iMax |
117 |
rmspp = rmspp |
118 |
& + (totPhiHyd(i,j,k,bi,bj)-oldPhi(i,j))**2 |
119 |
count = count + maskC(i,j,k,bi,bj) |
120 |
ENDDO |
121 |
ENDDO |
122 |
C -- end k loop |
123 |
ENDDO |
124 |
ENDDO |
125 |
ENDDO |
126 |
Cml WRITE(msgBuf,'(I10.10)') npiter |
127 |
Cml CALL WRITE_FLD_XYZ_RL( 'POLD.',msgBuf,pold,npiter,myThid) |
128 |
Cml CALL WRITE_FLD_XYZ_RL( 'PINI.',msgBuf,pressure,npiter,myThid) |
129 |
_GLOBAL_SUM_R8( rmspp, myThid ) |
130 |
_GLOBAL_SUM_R8( count, myThid ) |
131 |
IF ( count .EQ. 0. ) THEN |
132 |
rmspp = 0. _d 0 |
133 |
ELSE |
134 |
rmspp = sqrt(rmspp/count) |
135 |
ENDIF |
136 |
WRITE(msgBuf,'(a,i2,a,e20.13)') |
137 |
& 'Iteration ', npiter, ', RMS-difference = ', rmspp |
138 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
139 |
& SQUEEZE_RIGHT , 1) |
140 |
|
141 |
C -- end npiter loop |
142 |
ENDDO |
143 |
C print some diagnostics |
144 |
IF ( rmspp .ne. 0. ) THEN |
145 |
IF ( rmspp .EQ. rmsppold ) THEN |
146 |
WRITE(msgBuf,'(A)') |
147 |
& 'Initial hydrostatic pressure did not converge perfectly,' |
148 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
149 |
& SQUEEZE_RIGHT , 1) |
150 |
WRITE(msgBuf,'(A)') |
151 |
& 'but the RMS-difference is constant, indicating that the' |
152 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
153 |
& SQUEEZE_RIGHT , 1) |
154 |
WRITE(msgBuf,'(A)') |
155 |
& 'algorithm converged within machine precision.' |
156 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
157 |
& SQUEEZE_RIGHT , 1) |
158 |
ELSE |
159 |
WRITE(msgBuf,'(A,I2,A)') |
160 |
& 'Initial hydrostatic pressure did not converge after ', |
161 |
& npiter-1, ' steps' |
162 |
CALL PRINT_ERROR( msgBuf, myThid ) |
163 |
STOP 'ABNORMAL END: S/R INI_PRESSURE' |
164 |
ENDIF |
165 |
ENDIF |
166 |
WRITE(msgBuf,'(A)') |
167 |
& 'Initial hydrostatic pressure converged.' |
168 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
169 |
& SQUEEZE_RIGHT , 1) |
170 |
|
171 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
172 |
|
173 |
c-- else of if useDynP_inEos_Zc |
174 |
ELSE |
175 |
C print a message and DO nothing |
176 |
WRITE(msgBuf,'(A,A)') |
177 |
& 'Pressure is predetermined for buoyancyRelation ', |
178 |
& buoyancyRelation(1:11) |
179 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
180 |
& SQUEEZE_RIGHT , 1) |
181 |
|
182 |
ENDIF |
183 |
|
184 |
WRITE(msgBuf,'(A)') ' ' |
185 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
186 |
& SQUEEZE_RIGHT , 1) |
187 |
|
188 |
RETURN |
189 |
END |