/[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.9 - (show annotations) (download)
Fri Aug 4 00:45:47 2006 UTC (17 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58t_post, checkpoint58o_post, checkpoint58p_post, checkpoint58q_post, mitgcm_mapl_00, checkpoint58r_post, checkpoint58n_post, checkpoint58u_post, checkpoint58s_post
Changes since 1.8: +12 -1 lines
print only once (master thread)

1 C $Header: /u/gcmpack/MITgcm/model/src/ini_pressure.F,v 1.8 2006/07/25 22:22:08 jmc 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 _BEGIN_MASTER( myThid )
61 WRITE(msgBuf,'(a)')
62 & 'Start initial hydrostatic pressure computation'
63 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
64 & SQUEEZE_RIGHT , myThid)
65 _END_MASTER( myThid )
66
67 DO bj = myByLo(myThid), myByHi(myThid)
68 DO bi = myBxLo(myThid), myBxHi(myThid)
69 DO k = 1,Nr
70 DO j=1-OLy,sNy+OLy
71 DO i=1-OLx,sNx+OLx
72 totPhiHyd(i,j,k,bi,bj) = 0. _d 0
73 ENDDO
74 ENDDO
75 ENDDO
76 ENDDO
77 ENDDO
78
79 IF ( useDynP_inEos_Zc ) THEN
80
81 #ifndef ALLOW_AUTODIFF_TAMC
82 cph-- deal with this iterative loop for AD once it will
83 cph-- really be needed;
84 cph-- would need storing of totPhiHyd for each npiter
85
86 rmspp = 1. _d 0
87 rmsppold = 0. _d 0
88 npiter = 0
89
90 C iterate pressure p = integral of (g*rho(p)*dz)
91 DO npiter= 1, 15
92 rmsppold = rmspp
93 rmspp = 0. _d 0
94 count = 0.
95 DO bj = myByLo(myThid), myByHi(myThid)
96 DO bi = myBxLo(myThid), myBxHi(myThid)
97 DO j=1-OLy,sNy+OLy
98 DO i=1-OLx,sNx+OLx
99 phiHydF(i,j) = 0. _d 0
100 ENDDO
101 ENDDO
102 DO k = 1, Nr
103 C for each level save old pressure and compute new pressure
104 DO j=jMin,jMax
105 DO i=iMin,iMax
106 oldPhi(i,j) = totPhiHyd(i,j,k,bi,bj)
107 ENDDO
108 ENDDO
109 CALL CALC_PHI_HYD(
110 I bi, bj, iMin, iMax, jMin, jMax, k,
111 I theta, salt,
112 U phiHydF,
113 O phiHydC, dPhiHydX, dPhiHydY,
114 I startTime, nIter0, myThid)
115 C compute convergence criterion
116 DO j=jMin,jMax
117 DO i=iMin,iMax
118 rmspp = rmspp
119 & + (totPhiHyd(i,j,k,bi,bj)-oldPhi(i,j))**2
120 count = count + maskC(i,j,k,bi,bj)
121 ENDDO
122 ENDDO
123 C -- end k loop
124 ENDDO
125 ENDDO
126 ENDDO
127 Cml WRITE(msgBuf,'(I10.10)') npiter
128 Cml CALL WRITE_FLD_XYZ_RL( 'POLD.',msgBuf,pold,npiter,myThid)
129 Cml CALL WRITE_FLD_XYZ_RL( 'PINI.',msgBuf,pressure,npiter,myThid)
130 _GLOBAL_SUM_R8( rmspp, myThid )
131 _GLOBAL_SUM_R8( count, myThid )
132 IF ( count .EQ. 0. ) THEN
133 rmspp = 0. _d 0
134 ELSE
135 rmspp = sqrt(rmspp/count)
136 ENDIF
137 _BEGIN_MASTER( myThid )
138 WRITE(msgBuf,'(a,i2,a,e20.13)')
139 & 'Iteration ', npiter, ', RMS-difference = ', rmspp
140 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
141 & SQUEEZE_RIGHT , myThid)
142 _END_MASTER( myThid )
143
144 C -- end npiter loop
145 ENDDO
146 C print some diagnostics
147 _BEGIN_MASTER( myThid )
148 IF ( rmspp .ne. 0. ) THEN
149 IF ( rmspp .EQ. rmsppold ) THEN
150 WRITE(msgBuf,'(A)')
151 & 'Initial hydrostatic pressure did not converge perfectly,'
152 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
153 & SQUEEZE_RIGHT , myThid)
154 WRITE(msgBuf,'(A)')
155 & 'but the RMS-difference is constant, indicating that the'
156 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
157 & SQUEEZE_RIGHT , myThid)
158 WRITE(msgBuf,'(A)')
159 & 'algorithm converged within machine precision.'
160 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
161 & SQUEEZE_RIGHT , myThid)
162 ELSE
163 WRITE(msgBuf,'(A,I2,A)')
164 & 'Initial hydrostatic pressure did not converge after ',
165 & npiter-1, ' steps'
166 CALL PRINT_ERROR( msgBuf, myThid )
167 STOP 'ABNORMAL END: S/R INI_PRESSURE'
168 ENDIF
169 ENDIF
170 WRITE(msgBuf,'(A)')
171 & 'Initial hydrostatic pressure converged.'
172 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
173 & SQUEEZE_RIGHT , myThid)
174 _END_MASTER( myThid )
175
176 #endif /* ALLOW_AUTODIFF_TAMC */
177
178 c-- else of if useDynP_inEos_Zc
179 ELSE
180 C print a message and DO nothing
181 _BEGIN_MASTER( myThid )
182 WRITE(msgBuf,'(A,A)')
183 & 'Pressure is predetermined for buoyancyRelation ',
184 & buoyancyRelation(1:11)
185 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
186 & SQUEEZE_RIGHT , myThid)
187 _END_MASTER( myThid )
188
189 ENDIF
190
191 _BEGIN_MASTER( myThid )
192 WRITE(msgBuf,'(A)') ' '
193 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
194 & SQUEEZE_RIGHT , myThid)
195 _END_MASTER( myThid )
196
197 RETURN
198 END

  ViewVC Help
Powered by ViewVC 1.1.22