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

Annotation of /MITgcm/model/src/find_alpha.F

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


Revision 1.10 - (hide annotations) (download)
Wed Sep 18 16:30:34 2002 UTC (21 years, 8 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint46l_pre, checkpoint46j_pre, checkpoint46j_post, checkpoint46k_post, checkpoint46h_pre, checkpoint46g_post, checkpoint46i_post, checkpoint46h_post
Changes since 1.9: +9 -13 lines
cleaned up argument list

1 mlosch 1.10 C $Header: /u/gcmpack/MITgcm/model/src/find_alpha.F,v 1.9 2002/09/05 20:49:33 mlosch Exp $
2 cnh 1.5 C $Name: $
3 adcroft 1.1
4     #include "CPP_OPTIONS.h"
5     #define USE_FACTORIZED_POLY
6    
7 cnh 1.5 CBOP
8     C !ROUTINE: FIND_ALPHA
9     C !INTERFACE:
10     SUBROUTINE FIND_ALPHA (
11 mlosch 1.10 I bi, bj, iMin, iMax, jMin, jMax, k, kRef,
12 adcroft 1.1 O alphaloc )
13    
14 cnh 1.5 C !DESCRIPTION: \bv
15     C *==========================================================*
16     C | o SUBROUTINE FIND_ALPHA
17     C | Calculates [drho(S,T,z) / dT] of a horizontal slice
18     C *==========================================================*
19     C |
20     C | k - is the Theta/Salt level
21     C | kRef - determines pressure reference level
22     C | (not used in 'LINEAR' mode)
23     C |
24     C | alphaloc - drho / dT (kg/m^3/C)
25     C |
26     C *==========================================================*
27     C \ev
28    
29     C !USES:
30     IMPLICIT NONE
31 adcroft 1.1 c Common
32     #include "SIZE.h"
33     #include "DYNVARS.h"
34     #include "EEPARAMS.h"
35     #include "PARAMS.h"
36 mlosch 1.6 #include "EOS.h"
37     #include "GRID.h"
38 adcroft 1.1
39 cnh 1.5 C !INPUT/OUTPUT PARAMETERS:
40 adcroft 1.1 c Arguments
41     integer bi,bj,iMin,iMax,jMin,jMax
42     integer k ! Level of Theta/Salt slice
43     integer kRef ! Pressure reference level
44     _RL alphaloc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
45    
46 cnh 1.5 C !LOCAL VARIABLES:
47 adcroft 1.1 c Local
48     integer i,j
49     _RL refTemp,refSalt,tP,sP
50 mlosch 1.9 _RL t1, t2, t3, s1, s3o2, p1, p2, sp5, p1t1
51 mlosch 1.6 _RL drhoP0dtheta, drhoP0dthetaFresh, drhoP0dthetaSalt
52     _RL dKdtheta, dKdthetaFresh, dKdthetaSalt, dKdthetaPres
53     _RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
54     _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
55 mlosch 1.9 _RL dnum_dtheta, dden_dtheta
56     _RL rhoDen(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
57     _RL rhoLoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
58 mlosch 1.6 integer myThid
59 cnh 1.5 CEOP
60 adcroft 1.1
61 mlosch 1.6 Cml stop 'myThid is not properly defined in this subroutine'
62    
63 mlosch 1.9 if (equationOfState.eq.'LINEAR') then
64 adcroft 1.1
65     do j=jMin,jMax
66     do i=iMin,iMax
67     alphaloc(i,j) = -rhonil * tAlpha
68     enddo
69     enddo
70    
71 mlosch 1.9 elseif (equationOfState.eq.'POLY3') then
72 adcroft 1.1
73     refTemp=eosRefT(kRef)
74     refSalt=eosRefS(kRef)
75    
76     do j=jMin,jMax
77     do i=iMin,iMax
78     tP=theta(i,j,k,bi,bj)-refTemp
79     sP=salt(i,j,k,bi,bj)-refSalt
80     #ifdef USE_FACTORIZED_POLY
81     alphaloc(i,j) =
82     & ( eosC(6,kRef)
83     & *tP*3.
84     & +(eosC(7,kRef)*sP + eosC(3,kRef))*2.
85     & )*tP
86     & +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef)
87     &
88     #else
89     alphaloc(i,j) =
90     & eosC(1,kRef) +
91     & eosC(3,kRef)*tP*2. +
92     & eosC(4,kRef) *sP +
93     & eosC(6,kRef)*tP*tP*3. +
94     & eosC(7,kRef)*tP*2. *sP +
95     & eosC(8,kRef) *sP*sP
96     #endif
97     enddo
98     enddo
99    
100 mlosch 1.9 elseif ( equationOfState(1:5).eq.'JMD95'
101     & .or. equationOfState.eq.'UNESCO' ) then
102 mlosch 1.6 C nonlinear equation of state in pressure coordinates
103    
104     CALL FIND_RHOP0(
105     I bi, bj, iMin, iMax, jMin, jMax, k,
106     I theta, salt,
107     O rhoP0,
108     I myThid )
109    
110     CALL FIND_BULKMOD(
111 adcroft 1.7 I bi, bj, iMin, iMax, jMin, jMax, k, kRef,
112 mlosch 1.6 I theta, salt,
113     O bulkMod,
114     I myThid )
115    
116     do j=jMin,jMax
117     do i=iMin,iMax
118    
119     C abbreviations
120 mlosch 1.9 t1 = theta(i,j,k,bi,bj)
121     t2 = t1*t1
122     t3 = t2*t1
123 mlosch 1.6
124 mlosch 1.9 s1 = salt(i,j,k,bi,bj)
125     if ( s1 .lt. 0. _d 0 ) then
126 mlosch 1.6 C issue a warning
127 mlosch 1.10 write(*,'(a,i3,a,i3,a,i3,a,i3,a,i3,a,e13.5)')
128 mlosch 1.9 & ' FIND_ALPHA: WARNING, salinity(',
129     & i,',',j,',',k,',',bi,',',bj,') = ', s1
130     s1 = 0. _d 0
131 mlosch 1.6 end if
132 mlosch 1.9 s3o2 = sqrt(s1*s1*s1)
133 mlosch 1.6
134 mlosch 1.9 p1 = pressure(i,j,kRef,bi,bj)*SItoBar
135     p2 = p1*p1
136 mlosch 1.6
137     C d(rho)/d(theta)
138     C of fresh water at p = 0
139     drhoP0dthetaFresh =
140     & eosJMDCFw(2)
141 mlosch 1.9 & + 2.*eosJMDCFw(3)*t1
142 mlosch 1.6 & + 3.*eosJMDCFw(4)*t2
143     & + 4.*eosJMDCFw(5)*t3
144 mlosch 1.9 & + 5.*eosJMDCFw(6)*t3*t1
145 mlosch 1.6 C of salt water at p = 0
146     drhoP0dthetaSalt =
147 mlosch 1.9 & s1*(
148 mlosch 1.6 & eosJMDCSw(2)
149 mlosch 1.9 & + 2.*eosJMDCSw(3)*t1
150 mlosch 1.6 & + 3.*eosJMDCSw(4)*t2
151     & + 4.*eosJMDCSw(5)*t3
152     & )
153     & + s3o2*(
154     & + eosJMDCSw(7)
155 mlosch 1.9 & + 2.*eosJMDCSw(8)*t1
156 mlosch 1.6 & )
157     C d(bulk modulus)/d(theta)
158     C of fresh water at p = 0
159     dKdthetaFresh =
160     & eosJMDCKFw(2)
161 mlosch 1.9 & + 2.*eosJMDCKFw(3)*t1
162 mlosch 1.6 & + 3.*eosJMDCKFw(4)*t2
163     & + 4.*eosJMDCKFw(5)*t3
164     C of sea water at p = 0
165     dKdthetaSalt =
166 mlosch 1.9 & s1*( eosJMDCKSw(2)
167     & + 2.*eosJMDCKSw(3)*t1
168 mlosch 1.6 & + 3.*eosJMDCKSw(4)*t2
169     & )
170     & + s3o2*( eosJMDCKSw(6)
171 mlosch 1.9 & + 2.*eosJMDCKSw(7)*t1
172 mlosch 1.6 & )
173     C of sea water at p
174     dKdthetaPres =
175 mlosch 1.9 & p1*( eosJMDCKP(2)
176     & + 2.*eosJMDCKP(3)*t1
177 mlosch 1.6 & + 3.*eosJMDCKP(4)*t2
178     & )
179 mlosch 1.9 & + p1*s1*( eosJMDCKP(6)
180     & + 2.*eosJMDCKP(7)*t1
181 mlosch 1.6 & )
182     & + p2*( eosJMDCKP(10)
183 mlosch 1.9 & + 2.*eosJMDCKP(11)*t1
184 mlosch 1.6 & )
185 mlosch 1.9 & + p2*s1*( eosJMDCKP(13)
186     & + 2.*eosJMDCKP(14)*t1
187 mlosch 1.6 & )
188    
189     drhoP0dtheta = drhoP0dthetaFresh
190     & + drhoP0dthetaSalt
191     dKdtheta = dKdthetaFresh
192     & + dKdthetaSalt
193     & + dKdthetaPres
194     alphaloc(i,j) =
195     & ( bulkmod(i,j)**2*drhoP0dtheta
196 mlosch 1.9 & - bulkmod(i,j)*p1*drhoP0dtheta
197     & - rhoP0(i,j)*p1*dKdtheta )
198     & /( bulkmod(i,j) - p1 )**2
199 mlosch 1.6
200    
201     end do
202     end do
203 mlosch 1.9 elseif ( equationOfState.eq.'MDJWF' ) then
204    
205     CALL FIND_RHO( bi, bj, iMin, iMax, jMin, jMax, k, kRef,
206 mlosch 1.10 & theta, salt, rhoLoc, myThid )
207 mlosch 1.9 CALL FIND_RHODEN( bi, bj, iMin, iMax, jMin, jMax, k, kRef,
208     & theta, salt, rhoDen, myThid )
209    
210     do j=jMin,jMax
211     do i=iMin,iMax
212     t1 = theta(i,j,k,bi,bj)
213     t2 = t1*t1
214     s1 = salt(i,j,k,bi,bj)
215     if ( s1 .lt. 0. _d 0 ) then
216     C issue a warning
217 mlosch 1.10 write(*,'(a,i3,a,i3,a,i3,a,i3,a,i3,a,e13.5)')
218 mlosch 1.9 & ' FIND_ALPHA: WARNING, salinity(',
219     & i,',',j,',',k,',',bi,',',bj,') = ', s1
220     s1 = 0. _d 0
221     end if
222     sp5 = sqrt(s1)
223    
224     p1 = pressure(i,j,kRef,bi,bj)*SItodBar
225     p1t1 = p1*t1
226    
227     dnum_dtheta = eosMDJWFnum(1)
228     & + t1*(2.*eosMDJWFnum(2) + 3.*eosMDJWFnum(3)*t1)
229     & + eosMDJWFnum(5)*s1
230     & + p1t1*(2.*eosMDJWFnum(8) + 2.*eosMDJWFnum(11)*p1)
231    
232     dden_dtheta = eosMDJWFden(1)
233     & + t1*(2.*eosMDJWFden(2)
234     & + t1*(3.*eosMDJWFden(3)
235     & + 4.*eosMDJWFden(4)*t1 ) )
236     & + s1*(eosMDJWFden(6)
237     & + t1*(3.*eosMDJWFden(7)*t1
238     & + 2.*eosMDJWFden(9)*sp5 ) )
239     & + p1*p1*(3.*eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1)
240    
241     alphaLoc(i,j) = rhoDen(i,j)*(dnum_dtheta
242     & - rhoLoc(i,j)*dden_dtheta)
243    
244     end do
245     end do
246    
247 adcroft 1.1 else
248 mlosch 1.9 write(*,*) 'FIND_ALPHA: equationOfState = ',equationOfState
249     stop 'FIND_ALPHA: "equationOfState" has illegal value'
250 adcroft 1.1 endif
251    
252     return
253     end
254    
255     subroutine FIND_BETA (
256 mlosch 1.10 I bi, bj, iMin, iMax, jMin, jMax, k, kRef,
257 adcroft 1.1 O betaloc )
258     C /==========================================================\
259     C | o SUBROUTINE FIND_BETA |
260     C | Calculates [drho(S,T,z) / dS] of a horizontal slice |
261     C |==========================================================|
262     C | |
263     C | k - is the Theta/Salt level |
264     C | kRef - determines pressure reference level |
265     C | (not used in 'LINEAR' mode) |
266     C | |
267     C | betaloc - drho / dS (kg/m^3/PSU) |
268     C | |
269     C \==========================================================/
270     implicit none
271    
272     c Common
273     #include "SIZE.h"
274     #include "DYNVARS.h"
275     #include "EEPARAMS.h"
276     #include "PARAMS.h"
277 mlosch 1.6 #include "EOS.h"
278     #include "GRID.h"
279 adcroft 1.1
280     c Arguments
281     integer bi,bj,iMin,iMax,jMin,jMax
282     integer k ! Level of Theta/Salt slice
283     integer kRef ! Pressure reference level
284     _RL betaloc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
285    
286     c Local
287     integer i,j
288     _RL refTemp,refSalt,tP,sP
289 mlosch 1.9 _RL t1, t2, t3, s1, s3o2, p1, sp5, p1t1
290 mlosch 1.6 _RL drhoP0dS
291     _RL dKdS, dKdSSalt, dKdSPres
292     _RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
293     _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
294 mlosch 1.9 _RL dnum_dsalt, dden_dsalt
295     _RL rhoDen(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
296     _RL rhoLoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
297 mlosch 1.6 integer myThid
298     CEOP
299    
300     Cml stop 'myThid is not properly defined in this subroutine'
301 adcroft 1.1
302 mlosch 1.9 if (equationOfState.eq.'LINEAR') then
303 adcroft 1.1
304     do j=jMin,jMax
305     do i=iMin,iMax
306     betaloc(i,j) = rhonil * sBeta
307     enddo
308     enddo
309    
310 mlosch 1.9 elseif (equationOfState.eq.'POLY3') then
311 adcroft 1.1
312     refTemp=eosRefT(kRef)
313     refSalt=eosRefS(kRef)
314    
315     do j=jMin,jMax
316     do i=iMin,iMax
317     tP=theta(i,j,k,bi,bj)-refTemp
318     sP=salt(i,j,k,bi,bj)-refSalt
319     #ifdef USE_FACTORIZED_POLY
320     betaloc(i,j) =
321     & ( eosC(9,kRef)*sP*3. + eosC(5,kRef)*2. )*sP + eosC(2,kRef)
322     & + ( eosC(7,kRef)*tP
323     & +eosC(8,kRef)*sP*2. + eosC(4,kRef)
324     & )*tP
325     #else
326     betaloc(i,j) =
327     & eosC(2,kRef) +
328     & eosC(4,kRef)*tP +
329     & eosC(5,kRef) *sP*2. +
330     & eosC(7,kRef)*tP*tP +
331     & eosC(8,kRef)*tP *sP*2. +
332     & eosC(9,kRef) *sP*sP*3.
333     #endif
334     enddo
335     enddo
336    
337 mlosch 1.9 elseif ( equationOfState(1:5).eq.'JMD95'
338     & .or. equationOfState.eq.'UNESCO' ) then
339 mlosch 1.6 C nonlinear equation of state in pressure coordinates
340    
341     CALL FIND_RHOP0(
342     I bi, bj, iMin, iMax, jMin, jMax, k,
343     I theta, salt,
344     O rhoP0,
345     I myThid )
346    
347     CALL FIND_BULKMOD(
348 adcroft 1.7 I bi, bj, iMin, iMax, jMin, jMax, k, kRef,
349 mlosch 1.6 I theta, salt,
350     O bulkMod,
351     I myThid )
352    
353     do j=jMin,jMax
354     do i=iMin,iMax
355    
356     C abbreviations
357 mlosch 1.9 t1 = theta(i,j,k,bi,bj)
358     t2 = t1*t1
359     t3 = t2*t1
360 mlosch 1.6
361 mlosch 1.9 s1 = salt(i,j,k,bi,bj)
362     if ( s1 .lt. 0. _d 0 ) then
363 mlosch 1.6 C issue a warning
364 mlosch 1.10 write(*,'(a,i3,a,i3,a,i3,a,i3,a,i3,a,e13.5)')
365 mlosch 1.9 & ' FIND_BETA: WARNING, salinity(',
366     & i,',',j,',',k,',',bi,',',bj,') = ', s1
367     s1 = 0. _d 0
368 mlosch 1.6 end if
369 mlosch 1.9 s3o2 = 1.5*sqrt(s1)
370 mlosch 1.6
371 mlosch 1.9 p1 = pressure(i,j,kRef,bi,bj)*SItoBar
372 mlosch 1.6
373     C d(rho)/d(S)
374     C of fresh water at p = 0
375     drhoP0dS = 0. _d 0
376     C of salt water at p = 0
377     drhoP0dS = drhoP0dS
378     & + eosJMDCSw(1)
379 mlosch 1.9 & + eosJMDCSw(2)*t1
380 mlosch 1.6 & + eosJMDCSw(3)*t2
381     & + eosJMDCSw(4)*t3
382 mlosch 1.9 & + eosJMDCSw(5)*t3*t1
383 mlosch 1.6 & + s3o2*(
384     & eosJMDCSw(6)
385 mlosch 1.9 & + eosJMDCSw(7)*t1
386 mlosch 1.6 & + eosJMDCSw(8)*t2
387     & )
388 mlosch 1.9 & + 2*eosJMDCSw(9)*s1
389 mlosch 1.6 C d(bulk modulus)/d(S)
390     C of fresh water at p = 0
391     dKdS = 0. _d 0
392     C of sea water at p = 0
393     dKdSSalt =
394     & eosJMDCKSw(1)
395 mlosch 1.9 & + eosJMDCKSw(2)*t1
396 mlosch 1.6 & + eosJMDCKSw(3)*t2
397     & + eosJMDCKSw(4)*t3
398     & + s3o2*( eosJMDCKSw(5)
399 mlosch 1.9 & + eosJMDCKSw(6)*t1
400 mlosch 1.6 & + eosJMDCKSw(7)*t2
401     & )
402    
403     C of sea water at p
404     dKdSPres =
405 mlosch 1.9 & p1*( eosJMDCKP(5)
406     & + eosJMDCKP(6)*t1
407 mlosch 1.6 & + eosJMDCKP(7)*t2
408     & )
409 mlosch 1.9 & + s3o2*p1*eosJMDCKP(8)
410     & + p1*p1*( eosJMDCKP(12)
411     & + eosJMDCKP(13)*t1
412 mlosch 1.6 & + eosJMDCKP(14)*t2
413     & )
414    
415     dKdS = dKdSSalt + dKdSPres
416    
417     betaloc(i,j) =
418     & ( bulkmod(i,j)**2*drhoP0dS
419 mlosch 1.9 & - bulkmod(i,j)*p1*drhoP0dS
420     & - rhoP0(i,j)*p1*dKdS )
421     & /( bulkmod(i,j) - p1 )**2
422 mlosch 1.6
423    
424     end do
425     end do
426 mlosch 1.9 elseif ( equationOfState.eq.'MDJWF' ) then
427    
428     CALL FIND_RHO( bi, bj, iMin, iMax, jMin, jMax, k, kRef,
429 mlosch 1.10 & theta, salt, rhoLoc, myThid )
430 mlosch 1.9 CALL FIND_RHODEN( bi, bj, iMin, iMax, jMin, jMax, k, kRef,
431     & theta, salt, rhoDen, myThid )
432    
433     do j=jMin,jMax
434     do i=iMin,iMax
435     t1 = theta(i,j,k,bi,bj)
436     t2 = t1*t1
437     s1 = salt(i,j,k,bi,bj)
438     if ( s1 .lt. 0. _d 0 ) then
439     C issue a warning
440 mlosch 1.10 write(*,'(a,i3,a,i3,a,i3,a,i3,a,i3,a,e13.5)')
441 mlosch 1.9 & ' FIND_BETA: WARNING, salinity(',
442     & i,',',j,',',k,',',bi,',',bj,') = ', s1
443     s1 = 0. _d 0
444     end if
445     sp5 = sqrt(s1)
446    
447     p1 = pressure(i,j,k,bi,bj)*SItodBar
448     p1t1 = p1*t1
449    
450     dnum_dsalt = eosMDJWFnum(4)
451     & + eosMDJWFnum(5)*t1
452     & + 2.*eosMDJWFnum(6)*s1 + eosMDJWFnum(9)*p1
453     dden_dsalt = eosMDJWFden(5)
454     & + t1*( eosMDJWFden(6) + eosMDJWFden(7)*t2 )
455     & + 1.5*sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2)
456    
457     betaLoc(i,j) = rhoDen(i,j)*( dnum_dsalt
458     & - rhoLoc(i,j)*dden_dsalt )
459    
460     end do
461     end do
462    
463 adcroft 1.1 else
464 mlosch 1.9 write(*,*) 'FIND_BETA: equationOfState = ',equationOfState
465     stop 'FIND_BETA: "equationOfState" has illegal value'
466 adcroft 1.1 endif
467    
468     return
469     end

  ViewVC Help
Powered by ViewVC 1.1.22