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

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

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


Revision 1.4 - (show annotations) (download)
Tue Feb 18 15:34:25 2003 UTC (21 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint48i_post, checkpoint48h_post, checkpoint49, checkpoint48g_post
Branch point for: ecco-branch
Changes since 1.3: +34 -29 lines
o compute locally the pressure for use in EOS : UNESCO, JMD95P or MDJWF
o store total Potential in totPhyHyd for diagnostic (DYNVARS.h)
o fix restart and overlap Pb when using Z-coord and EOS funct. of P

1 C $Header:
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 rmspp = 1. _d 0
81 rmsppold = 0. _d 0
82 npiter = 0
83
84 C iterate pressure p = integral of (g*rho(p)*dz)
85 DO npiter= 1, 15
86 rmsppold = rmspp
87 rmspp = 0. _d 0
88 count = 0.
89 DO bj = myByLo(myThid), myByHi(myThid)
90 DO bi = myBxLo(myThid), myBxHi(myThid)
91 DO j=1-OLy,sNy+OLy
92 DO i=1-OLx,sNx+OLx
93 phiHydF(i,j) = 0. _d 0
94 ENDDO
95 ENDDO
96 DO k = 1, Nr
97 C for each level save old pressure and compute new pressure
98 DO j=jMin,jMax
99 DO i=iMin,iMax
100 oldPhi(i,j) = totPhiHyd(i,j,k,bi,bj)
101 ENDDO
102 ENDDO
103 CALL CALC_PHI_HYD(
104 I bi, bj, iMin, iMax, jMin, jMax, k,
105 I theta, salt,
106 U phiHydF,
107 O phiHydC, dPhiHydX, dPhiHydY,
108 I startTime, nIter0, myThid)
109 C compute convergence criterion
110 DO j=jMin,jMax
111 DO i=iMin,iMax
112 rmspp = rmspp
113 & + (totPhiHyd(i,j,k,bi,bj)-oldPhi(i,j))**2
114 count = count + maskC(i,j,k,bi,bj)
115 ENDDO
116 ENDDO
117 C -- end k loop
118 ENDDO
119 ENDDO
120 ENDDO
121 Cml WRITE(msgBuf,'(I10.10)') npiter
122 Cml CALL WRITE_FLD_XYZ_RL( 'POLD.',msgBuf,pold,npiter,myThid)
123 Cml CALL WRITE_FLD_XYZ_RL( 'PINI.',msgBuf,pressure,npiter,myThid)
124 _GLOBAL_SUM_R8( rmspp, myThid )
125 _GLOBAL_SUM_R8( count, myThid )
126 IF ( count .EQ. 0. ) THEN
127 rmspp = 0. _d 0
128 ELSE
129 rmspp = sqrt(rmspp/count)
130 ENDIF
131 WRITE(msgBuf,'(a,i2,a,e20.13)')
132 & 'Iteration ', npiter, ', RMS-difference = ', rmspp
133 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
134 & SQUEEZE_RIGHT , 1)
135
136 C -- end npiter loop
137 ENDDO
138 C print some diagnostics
139 IF ( rmspp .ne. 0. ) THEN
140 IF ( rmspp .EQ. rmsppold ) THEN
141 WRITE(msgBuf,'(A)')
142 & 'Initial hydrostatic pressure did not converge perfectly,'
143 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
144 & SQUEEZE_RIGHT , 1)
145 WRITE(msgBuf,'(A)')
146 & 'but the RMS-difference is constant, indicating that the'
147 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
148 & SQUEEZE_RIGHT , 1)
149 WRITE(msgBuf,'(A)')
150 & 'algorithm converged within machine precision.'
151 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
152 & SQUEEZE_RIGHT , 1)
153 ELSE
154 WRITE(msgBuf,'(A,I2,A)')
155 & 'Initial hydrostatic pressure did not converge after ',
156 & npiter-1, ' steps'
157 CALL PRINT_ERROR( msgBuf, myThid )
158 STOP 'ABNORMAL END: S/R INI_PRESSURE'
159 ENDIF
160 ENDIF
161 WRITE(msgBuf,'(A)')
162 & 'Initial hydrostatic pressure converged.'
163 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
164 & SQUEEZE_RIGHT , 1)
165
166 ELSE
167 C print a message and DO nothing
168 WRITE(msgBuf,'(A,A)')
169 & 'Pressure is predetermined for buoyancyRelation ',
170 & buoyancyRelation(1:11)
171 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
172 & SQUEEZE_RIGHT , 1)
173
174 ENDIF
175
176 WRITE(msgBuf,'(A)') ' '
177 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
178 & SQUEEZE_RIGHT , 1)
179
180 RETURN
181 END

  ViewVC Help
Powered by ViewVC 1.1.22