/[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.7 - (show annotations) (download)
Wed Apr 6 18:29:53 2005 UTC (19 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57g_pre, checkpoint57s_post, checkpoint58b_post, checkpoint57g_post, checkpoint57y_post, checkpoint57r_post, checkpoint57i_post, checkpoint58, checkpoint58f_post, checkpoint57n_post, checkpoint58d_post, checkpoint58a_post, checkpoint57z_post, checkpoint58m_post, checkpoint57l_post, checkpoint57t_post, checkpoint57v_post, checkpoint57f_post, checkpoint57h_pre, checkpoint57h_post, checkpoint57y_pre, checkpoint58e_post, checkpoint57p_post, checkpint57u_post, checkpoint57q_post, checkpoint58k_post, checkpoint58l_post, checkpoint57h_done, checkpoint57j_post, checkpoint58g_post, checkpoint58h_post, checkpoint58j_post, checkpoint57o_post, checkpoint57k_post, checkpoint57w_post, checkpoint58i_post, checkpoint57x_post, checkpoint58c_post
Changes since 1.6: +1 -3 lines
use baseTime as time origin ; DIFF_BASE_MULTIPLE replaces DIFFERENT_MULTIPLE

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

  ViewVC Help
Powered by ViewVC 1.1.22