/[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.3 - (show annotations) (download)
Sat Feb 8 02:09:20 2003 UTC (21 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint48e_post, checkpoint48d_pre, checkpoint48d_post, checkpoint48f_post
Changes since 1.2: +68 -60 lines
in preparation for r*:
 new S/R (calc_grad_phi_hyd.F) to compute Hydrostatic potential gradient.
 pass the 2 comp. of the grad. as arguments to momentum S/R.
for the moment, only used if it does not change the results.

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 ncount, iMin, iMax, jMin, jMax, npiter
46 _RL phiHyd(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
47 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
48 _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
49 _RL pold(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
50 _RL rmspp, rmsppold
51
52 CHARACTER*(MAX_LEN_MBUF) msgBuf
53 CEOP
54 iMin = 1-OLx
55 iMax = sNx+OLx
56 jMin = 1-OLy
57 jMax = sNy+OLy
58
59 WRITE(msgBuf,'(a)')
60 & 'Start initial hydrostatic pressure computation'
61 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
62 & SQUEEZE_RIGHT , 1)
63 C initialise pressure to water column thickness times a constant
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 phiHyd(i,j,k) = 0. _d 0
70 dPhiHydX(i,j) = 0. _d 0
71 dPhiHydY(i,j) = 0. _d 0
72 ENDDO
73 ENDDO
74 CALL STORE_PRESSURE( bi, bj, k, phiHyd, myThid )
75 ENDDO
76 ENDDO
77 ENDDO
78
79 IF ( buoyancyRelation .EQ. 'OCEANIC' ) THEN
80
81 rmspp = 1. _d 0
82 rmsppold = 0. _d 0
83 npiter = 0
84
85 C iterate pressure p = integral of (g*rho(p)*dz)
86 DO npiter = 1, 15
87 rmsppold = rmspp
88 rmspp = 0. _d 0
89 ncount = 0
90 DO bj = myByLo(myThid), myByHi(myThid)
91 DO bi = myBxLo(myThid), myBxHi(myThid)
92 DO k = 1, Nr
93 C for each level save old pressure and compute new pressure
94 DO j=jMin,jMax
95 DO i=iMin,iMax
96 pold(i,j,k) = pressure(i,j,k,bi,bj)
97 ENDDO
98 ENDDO
99 CALL CALC_PHI_HYD(
100 I bi, bj, iMin, iMax, jMin, jMax, k,
101 I theta, salt,
102 U phiHyd,
103 O dPhiHydX, dPhiHydY,
104 I startTime, nIter0, myThid)
105 CALL STORE_PRESSURE( bi, bj, k, phiHyd, myThid )
106 ENDDO
107 C compute convergence criterion
108 DO k = 1, Nr
109 DO j=jMin,jMax
110 DO i=iMin,iMax
111 rmspp = rmspp
112 & + (pressure(i,j,k,bi,bj)-pold(i,j,k))**2
113 ncount = ncount + int(maskC(i,j,k,bi,bj))
114 ENDDO
115 ENDDO
116 ENDDO
117 ENDDO
118 ENDDO
119 Cml WRITE(msgBuf,'(I10.10)') npiter
120 Cml CALL WRITE_FLD_XYZ_RL( 'POLD.',msgBuf,pold,npiter,myThid)
121 Cml CALL WRITE_FLD_XYZ_RL( 'PINI.',msgBuf,pressure,npiter,myThid)
122 IF ( ncount .EQ. 0 ) THEN
123 rmspp = 0. _d 0
124 ELSE
125 rmspp = sqrt(rmspp/ncount)
126 ENDIF
127 WRITE(msgBuf,'(a,i2,a,e20.13)')
128 & 'Iteration ', npiter, ', RMS-difference = ', rmspp
129 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
130 & SQUEEZE_RIGHT , 1)
131
132 ENDDO
133 C print some diagnostics
134 IF ( rmspp .ne. 0. ) THEN
135 IF ( rmspp .EQ. rmsppold ) THEN
136 WRITE(msgBuf,'(a)')
137 & 'Initial hydrostatic pressure did not converge perfectly,'
138 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
139 & SQUEEZE_RIGHT , 1)
140 WRITE(msgBuf,'(a)')
141 & 'but the RMS-difference is constant, indicating that the'
142 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
143 & SQUEEZE_RIGHT , 1)
144 WRITE(msgBuf,'(a)')
145 & 'algorithm converged within machine precision.'
146 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
147 & SQUEEZE_RIGHT , 1)
148 ELSE
149 WRITE(msgBuf,'(a,i2,a)')
150 & 'Initial hydrostatic pressure did not converge after ',
151 & npiter-1, ' steps'
152 CALL PRINT_ERROR( msgBuf, myThid )
153 STOP 'ABNORMAL END: S/R INI_PRESSURE'
154 ENDIF
155 ENDIF
156 WRITE(msgBuf,'(a)')
157 & 'Initial hydrostatic pressure converged.'
158 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
159 & SQUEEZE_RIGHT , 1)
160
161 ELSE
162 C print a message and DO nothing
163 WRITE(msgBuf,'(a,a)')
164 & 'Pressure is predetermined for buoyancyRelation ',
165 & buoyancyRelation(1:11)
166 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
167 & SQUEEZE_RIGHT , 1)
168
169 ENDIF
170
171 WRITE(msgBuf,'(a)') ' '
172 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
173 & SQUEEZE_RIGHT , 1)
174
175 RETURN
176 END

  ViewVC Help
Powered by ViewVC 1.1.22