/[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.18 - (hide annotations) (download)
Tue Jul 19 12:53:24 2011 UTC (12 years, 10 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint63g, checkpoint63h, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c
Changes since 1.17: +95 -1 lines
add code for TEOS-10 (www.teos-10.org, McDougall et al. 2011). Use
this eos with eosType = 'TEOS10', in data (PARM01). This eos implies
that THETA and SALT are "conservative temperature" and "absolute
salinity"

1 mlosch 1.18 C $Header: /u/gcmpack/MITgcm/model/src/find_alpha.F,v 1.17 2008/08/09 01:02:28 jmc 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 jmc 1.17 I bi, bj, iMin, iMax, jMin, jMax, k, kRef,
12     O alphaLoc,
13 mlosch 1.16 I myThid )
14 adcroft 1.1
15 cnh 1.5 C !DESCRIPTION: \bv
16     C *==========================================================*
17 jmc 1.17 C | o SUBROUTINE FIND_ALPHA
18     C | Calculates [drho(S,T,z) / dT] of a horizontal slice
19 cnh 1.5 C *==========================================================*
20 jmc 1.17 C |
21     C | k - is the Theta/Salt level
22     C | kRef - determines pressure reference level
23     C | (not used in 'LINEAR' mode)
24     C |
25     C | alphaLoc - drho / dT (kg/m^3/C)
26     C |
27 cnh 1.5 C *==========================================================*
28     C \ev
29    
30     C !USES:
31     IMPLICIT NONE
32 jmc 1.12 C === Global variables ===
33 adcroft 1.1 #include "SIZE.h"
34     #include "DYNVARS.h"
35     #include "EEPARAMS.h"
36     #include "PARAMS.h"
37 mlosch 1.6 #include "EOS.h"
38     #include "GRID.h"
39 adcroft 1.1
40 cnh 1.5 C !INPUT/OUTPUT PARAMETERS:
41 jmc 1.12 C == Routine arguments ==
42 mlosch 1.16 C k :: Level of Theta/Salt slice
43     C kRef :: Pressure reference level
44     c myThid :: thread number for this instance of the routine
45     INTEGER myThid
46 jmc 1.12 INTEGER bi,bj,iMin,iMax,jMin,jMax
47     INTEGER k
48     INTEGER kRef
49 jmc 1.17 _RL alphaLoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
50 adcroft 1.1
51 cnh 1.5 C !LOCAL VARIABLES:
52 jmc 1.17 C == Local variables ==
53 jmc 1.12 INTEGER i,j
54 adcroft 1.1 _RL refTemp,refSalt,tP,sP
55 mlosch 1.9 _RL t1, t2, t3, s1, s3o2, p1, p2, sp5, p1t1
56 mlosch 1.18 _RL ct, sa, sqrtsa, p
57 mlosch 1.6 _RL drhoP0dtheta, drhoP0dthetaFresh, drhoP0dthetaSalt
58     _RL dKdtheta, dKdthetaFresh, dKdthetaSalt, dKdthetaPres
59 jmc 1.12 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
60     _RL rhoP0 (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
61 mlosch 1.6 _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
62 mlosch 1.9 _RL dnum_dtheta, dden_dtheta
63 jmc 1.12 _RL rhoDen (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
64     _RL rhoLoc (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
65 cnh 1.5 CEOP
66 adcroft 1.1
67 mlosch 1.11 #ifdef CHECK_SALINITY_FOR_NEGATIVE_VALUES
68 jmc 1.17 c CALL LOOK_FOR_NEG_SALINITY(
69     c I iMin, iMax, jMin, jMax,
70     c U sFld,
71     c I k, bi, bj, myThid )
72 mlosch 1.11 #endif
73    
74 jmc 1.12 IF (equationOfState.EQ.'LINEAR') THEN
75 adcroft 1.1
76 jmc 1.12 DO j=jMin,jMax
77     DO i=iMin,iMax
78 jmc 1.17 alphaLoc(i,j) = -rhonil * tAlpha
79 jmc 1.12 ENDDO
80     ENDDO
81 jmc 1.17
82 jmc 1.12 ELSEIF (equationOfState.EQ.'POLY3') THEN
83 adcroft 1.1
84     refTemp=eosRefT(kRef)
85     refSalt=eosRefS(kRef)
86    
87 jmc 1.12 DO j=jMin,jMax
88     DO i=iMin,iMax
89 adcroft 1.1 tP=theta(i,j,k,bi,bj)-refTemp
90     sP=salt(i,j,k,bi,bj)-refSalt
91     #ifdef USE_FACTORIZED_POLY
92 jmc 1.17 alphaLoc(i,j) =
93 adcroft 1.1 & ( eosC(6,kRef)
94     & *tP*3.
95     & +(eosC(7,kRef)*sP + eosC(3,kRef))*2.
96     & )*tP
97     & +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef)
98 jmc 1.17 &
99 adcroft 1.1 #else
100 jmc 1.17 alphaLoc(i,j) =
101 adcroft 1.1 & eosC(1,kRef) +
102     & eosC(3,kRef)*tP*2. +
103     & eosC(4,kRef) *sP +
104     & eosC(6,kRef)*tP*tP*3. +
105     & eosC(7,kRef)*tP*2. *sP +
106     & eosC(8,kRef) *sP*sP
107     #endif
108 jmc 1.12 ENDDO
109     ENDDO
110 jmc 1.17
111 jmc 1.12 ELSEIF ( equationOfState(1:5).EQ.'JMD95'
112     & .OR. equationOfState.EQ.'UNESCO' ) THEN
113 mlosch 1.6 C nonlinear equation of state in pressure coordinates
114    
115 jmc 1.12 CALL PRESSURE_FOR_EOS(
116     I bi, bj, iMin, iMax, jMin, jMax, kRef,
117     O locPres,
118 jmc 1.17 I myThid )
119 jmc 1.12
120 mlosch 1.6 CALL FIND_RHOP0(
121 jmc 1.17 I iMin, iMax, jMin, jMax,
122     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
123 mlosch 1.6 O rhoP0,
124     I myThid )
125 jmc 1.17
126 mlosch 1.6 CALL FIND_BULKMOD(
127 jmc 1.17 I iMin, iMax, jMin, jMax, locPres,
128     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
129 mlosch 1.6 O bulkMod,
130     I myThid )
131    
132 jmc 1.12 DO j=jMin,jMax
133     DO i=iMin,iMax
134 mlosch 1.6
135     C abbreviations
136 jmc 1.17 t1 = theta(i,j,k,bi,bj)
137 mlosch 1.9 t2 = t1*t1
138     t3 = t2*t1
139 jmc 1.17
140 mlosch 1.9 s1 = salt(i,j,k,bi,bj)
141 jmc 1.14 IF ( s1 .GT. 0. _d 0 ) THEN
142 jmc 1.12 s3o2 = SQRT(s1*s1*s1)
143 jmc 1.14 ELSE
144     s1 = 0. _d 0
145     s3o2 = 0. _d 0
146     ENDIF
147 jmc 1.17
148 jmc 1.12 p1 = locPres(i,j)*SItoBar
149 mlosch 1.9 p2 = p1*p1
150 mlosch 1.6
151     C d(rho)/d(theta)
152     C of fresh water at p = 0
153 jmc 1.17 drhoP0dthetaFresh =
154 mlosch 1.6 & eosJMDCFw(2)
155 mlosch 1.9 & + 2.*eosJMDCFw(3)*t1
156 mlosch 1.6 & + 3.*eosJMDCFw(4)*t2
157     & + 4.*eosJMDCFw(5)*t3
158 mlosch 1.9 & + 5.*eosJMDCFw(6)*t3*t1
159 mlosch 1.6 C of salt water at p = 0
160 jmc 1.17 drhoP0dthetaSalt =
161 mlosch 1.9 & s1*(
162 mlosch 1.6 & eosJMDCSw(2)
163 mlosch 1.9 & + 2.*eosJMDCSw(3)*t1
164 mlosch 1.6 & + 3.*eosJMDCSw(4)*t2
165     & + 4.*eosJMDCSw(5)*t3
166     & )
167     & + s3o2*(
168     & + eosJMDCSw(7)
169 mlosch 1.9 & + 2.*eosJMDCSw(8)*t1
170 mlosch 1.6 & )
171     C d(bulk modulus)/d(theta)
172     C of fresh water at p = 0
173 jmc 1.17 dKdthetaFresh =
174 mlosch 1.6 & eosJMDCKFw(2)
175 mlosch 1.9 & + 2.*eosJMDCKFw(3)*t1
176 mlosch 1.6 & + 3.*eosJMDCKFw(4)*t2
177     & + 4.*eosJMDCKFw(5)*t3
178     C of sea water at p = 0
179     dKdthetaSalt =
180 mlosch 1.9 & s1*( eosJMDCKSw(2)
181     & + 2.*eosJMDCKSw(3)*t1
182 mlosch 1.6 & + 3.*eosJMDCKSw(4)*t2
183     & )
184     & + s3o2*( eosJMDCKSw(6)
185 mlosch 1.9 & + 2.*eosJMDCKSw(7)*t1
186 mlosch 1.6 & )
187     C of sea water at p
188     dKdthetaPres =
189 mlosch 1.9 & p1*( eosJMDCKP(2)
190     & + 2.*eosJMDCKP(3)*t1
191 mlosch 1.6 & + 3.*eosJMDCKP(4)*t2
192     & )
193 mlosch 1.9 & + p1*s1*( eosJMDCKP(6)
194     & + 2.*eosJMDCKP(7)*t1
195 mlosch 1.6 & )
196     & + p2*( eosJMDCKP(10)
197 mlosch 1.9 & + 2.*eosJMDCKP(11)*t1
198 mlosch 1.6 & )
199 mlosch 1.9 & + p2*s1*( eosJMDCKP(13)
200     & + 2.*eosJMDCKP(14)*t1
201 mlosch 1.6 & )
202    
203 jmc 1.17 drhoP0dtheta = drhoP0dthetaFresh
204 mlosch 1.6 & + drhoP0dthetaSalt
205 jmc 1.17 dKdtheta = dKdthetaFresh
206     & + dKdthetaSalt
207 mlosch 1.6 & + dKdthetaPres
208 jmc 1.17 alphaLoc(i,j) =
209 mlosch 1.6 & ( bulkmod(i,j)**2*drhoP0dtheta
210 mlosch 1.9 & - bulkmod(i,j)*p1*drhoP0dtheta
211     & - rhoP0(i,j)*p1*dKdtheta )
212     & /( bulkmod(i,j) - p1 )**2
213 jmc 1.17
214    
215 jmc 1.12 ENDDO
216     ENDDO
217     ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN
218    
219     CALL PRESSURE_FOR_EOS(
220     I bi, bj, iMin, iMax, jMin, jMax, kRef,
221     O locPres,
222 jmc 1.17 I myThid )
223 jmc 1.12
224 jmc 1.17 CALL FIND_RHONUM(
225     I iMin, iMax, jMin, jMax, locPres,
226     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
227     O rhoLoc,
228     I myThid )
229 heimbach 1.15
230 jmc 1.17 CALL FIND_RHODEN(
231     I iMin, iMax, jMin, jMax, locPres,
232     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
233     O rhoDen,
234     I myThid )
235 mlosch 1.9
236 jmc 1.12 DO j=jMin,jMax
237     DO i=iMin,iMax
238 mlosch 1.9 t1 = theta(i,j,k,bi,bj)
239     t2 = t1*t1
240     s1 = salt(i,j,k,bi,bj)
241 jmc 1.14 IF ( s1 .GT. 0. _d 0 ) THEN
242 jmc 1.17 sp5 = SQRT(s1)
243 jmc 1.14 ELSE
244     s1 = 0. _d 0
245     sp5 = 0. _d 0
246     ENDIF
247 mlosch 1.9
248 jmc 1.12 p1 = locPres(i,j)*SItodBar
249 mlosch 1.9 p1t1 = p1*t1
250 jmc 1.17
251 mlosch 1.9 dnum_dtheta = eosMDJWFnum(1)
252 jmc 1.17 & + t1*(2.*eosMDJWFnum(2) + 3.*eosMDJWFnum(3)*t1)
253     & + eosMDJWFnum(5)*s1
254 jmc 1.14 & + p1t1*(2.*eosMDJWFnum(8) + 2.*eosMDJWFnum(11)*p1)
255 jmc 1.17
256     dden_dtheta = eosMDJWFden(1)
257     & + t1*(2.*eosMDJWFden(2)
258     & + t1*(3.*eosMDJWFden(3)
259 mlosch 1.9 & + 4.*eosMDJWFden(4)*t1 ) )
260 jmc 1.17 & + s1*(eosMDJWFden(6)
261     & + t1*(3.*eosMDJWFden(7)*t1
262 mlosch 1.9 & + 2.*eosMDJWFden(9)*sp5 ) )
263     & + p1*p1*(3.*eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1)
264 jmc 1.17
265     alphaLoc(i,j) = rhoDen(i,j)*(dnum_dtheta
266 heimbach 1.15 & - (rhoLoc(i,j)*rhoDen(i,j))*dden_dtheta)
267 jmc 1.17
268 jmc 1.12 ENDDO
269     ENDDO
270 mlosch 1.9
271 mlosch 1.18 ELSEIF ( equationOfState.EQ.'TEOS10' ) THEN
272    
273     CALL PRESSURE_FOR_EOS(
274     I bi, bj, iMin, iMax, jMin, jMax, kRef,
275     O locPres,
276     I myThid )
277    
278     CALL FIND_RHOTEOS(
279     I iMin, iMax, jMin, jMax, locPres,
280     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
281     O rhoLoc, rhoDen,
282     I myThid )
283    
284     DO j=jMin,jMax
285     DO i=iMin,iMax
286     ct = theta(i,j,k,bi,bj)
287     sa = salt(i,j,k,bi,bj)
288     IF ( sa .GT. 0. _d 0 ) THEN
289     sqrtsa = SQRT(sa)
290     ELSE
291     sa = 0. _d 0
292     sqrtsa = 0. _d 0
293     ENDIF
294     p = locPres(i,j)*SItodBar
295    
296     dnum_dtheta = teos(02)
297     & + ct*(2.*teos(03) + 3.*teos(04)*ct)
298     & + sa*(teos(06) + 2.*teos(07)*ct
299     & + sqrtsa*(teos(09) + ct*(2.*teos(10) + 3.*teos(11)*ct)))
300     & + p*( teos(13) + 2.*teos(14)*ct + sa*2.*teos(16)
301     & + p*(teos(18) + 2.*teos(19)*ct))
302    
303     dden_dtheta = teos(22)
304     & + ct*(2.*teos(23) + ct*(3.*teos(24) + 4.*teos(25)*ct))
305     & + sa*(teos(27)
306     & + ct*(2.*teos(28) + ct*(3.*teos(29) + 4.*teos(30)*ct))
307     & + sqrtsa*(teos(32)
308     & + ct*(2.*teos(33) + ct*(3.*teos(34) + 4.*teos(35)*ct))))
309     & + p*(teos(38) + ct*(2.*teos(39) + 3.*teos(40)*ct)
310     & + teos(42)
311     & + p*(teos(44) + 2.*teos(45)*ct + teos(46)*sa
312     & + p*teos(48) ))
313    
314     alphaLoc(i,j) = rhoDen(i,j)*(dnum_dtheta
315     & - (rhoLoc(i,j)*rhoDen(i,j))*dden_dtheta)
316    
317     ENDDO
318     ENDDO
319    
320 jmc 1.12 ELSE
321     WRITE(*,*) 'FIND_ALPHA: equationOfState = ',equationOfState
322     STOP 'FIND_ALPHA: "equationOfState" has illegal value'
323     ENDIF
324 adcroft 1.1
325 jmc 1.17 RETURN
326 jmc 1.12 END
327 adcroft 1.1
328 jmc 1.12 SUBROUTINE FIND_BETA (
329 mlosch 1.10 I bi, bj, iMin, iMax, jMin, jMax, k, kRef,
330 jmc 1.17 O betaLoc,
331 mlosch 1.16 I myThid )
332 adcroft 1.1 C /==========================================================\
333     C | o SUBROUTINE FIND_BETA |
334     C | Calculates [drho(S,T,z) / dS] of a horizontal slice |
335     C |==========================================================|
336     C | |
337     C | k - is the Theta/Salt level |
338     C | kRef - determines pressure reference level |
339     C | (not used in 'LINEAR' mode) |
340     C | |
341 jmc 1.17 C | betaLoc - drho / dS (kg/m^3/PSU) |
342 adcroft 1.1 C | |
343     C \==========================================================/
344 jmc 1.12 IMPLICIT NONE
345 adcroft 1.1
346 jmc 1.12 C === Global variables ===
347 adcroft 1.1 #include "SIZE.h"
348     #include "DYNVARS.h"
349     #include "EEPARAMS.h"
350     #include "PARAMS.h"
351 mlosch 1.6 #include "EOS.h"
352     #include "GRID.h"
353 adcroft 1.1
354 jmc 1.12 C == Routine arguments ==
355 mlosch 1.16 C k :: Level of Theta/Salt slice
356     C kRef :: Pressure reference level
357     c myThid :: thread number for this instance of the routine
358     INTEGER myThid
359 jmc 1.12 INTEGER bi,bj,iMin,iMax,jMin,jMax
360     INTEGER k
361     INTEGER kRef
362 jmc 1.17 _RL betaLoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
363 adcroft 1.1
364 jmc 1.12 C == Local variables ==
365     INTEGER i,j
366 adcroft 1.1 _RL refTemp,refSalt,tP,sP
367 mlosch 1.9 _RL t1, t2, t3, s1, s3o2, p1, sp5, p1t1
368 mlosch 1.18 _RL ct, sa, sqrtsa, p
369 mlosch 1.6 _RL drhoP0dS
370     _RL dKdS, dKdSSalt, dKdSPres
371 jmc 1.12 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
372     _RL rhoP0 (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
373 mlosch 1.6 _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
374 mlosch 1.9 _RL dnum_dsalt, dden_dsalt
375 jmc 1.12 _RL rhoDen (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
376     _RL rhoLoc (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
377 mlosch 1.6 CEOP
378    
379 mlosch 1.11 #ifdef CHECK_SALINITY_FOR_NEGATIVE_VALUES
380 jmc 1.17 c CALL LOOK_FOR_NEG_SALINITY(
381     c I iMin, iMax, jMin, jMax,
382     c U sFld,
383     c I k, bi, bj, myThid )
384 mlosch 1.11 #endif
385    
386 jmc 1.12 IF (equationOfState.EQ.'LINEAR') THEN
387 adcroft 1.1
388 jmc 1.12 DO j=jMin,jMax
389     DO i=iMin,iMax
390 jmc 1.17 betaLoc(i,j) = rhonil * sBeta
391 jmc 1.12 ENDDO
392     ENDDO
393 jmc 1.17
394 jmc 1.12 ELSEIF (equationOfState.EQ.'POLY3') THEN
395 adcroft 1.1
396     refTemp=eosRefT(kRef)
397     refSalt=eosRefS(kRef)
398    
399 jmc 1.12 DO j=jMin,jMax
400     DO i=iMin,iMax
401 adcroft 1.1 tP=theta(i,j,k,bi,bj)-refTemp
402     sP=salt(i,j,k,bi,bj)-refSalt
403     #ifdef USE_FACTORIZED_POLY
404 jmc 1.17 betaLoc(i,j) =
405 adcroft 1.1 & ( eosC(9,kRef)*sP*3. + eosC(5,kRef)*2. )*sP + eosC(2,kRef)
406     & + ( eosC(7,kRef)*tP
407     & +eosC(8,kRef)*sP*2. + eosC(4,kRef)
408     & )*tP
409     #else
410 jmc 1.17 betaLoc(i,j) =
411 adcroft 1.1 & eosC(2,kRef) +
412     & eosC(4,kRef)*tP +
413     & eosC(5,kRef) *sP*2. +
414     & eosC(7,kRef)*tP*tP +
415     & eosC(8,kRef)*tP *sP*2. +
416     & eosC(9,kRef) *sP*sP*3.
417     #endif
418 jmc 1.12 ENDDO
419     ENDDO
420 adcroft 1.1
421 jmc 1.12 ELSEIF ( equationOfState(1:5).EQ.'JMD95'
422     & .OR. equationOfState.EQ.'UNESCO' ) THEN
423 mlosch 1.6 C nonlinear equation of state in pressure coordinates
424    
425 jmc 1.12 CALL PRESSURE_FOR_EOS(
426     I bi, bj, iMin, iMax, jMin, jMax, kRef,
427     O locPres,
428 jmc 1.17 I myThid )
429 jmc 1.12
430 mlosch 1.6 CALL FIND_RHOP0(
431 jmc 1.17 I iMin, iMax, jMin, jMax,
432     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
433 mlosch 1.6 O rhoP0,
434     I myThid )
435 jmc 1.17
436 mlosch 1.6 CALL FIND_BULKMOD(
437 jmc 1.17 I iMin, iMax, jMin, jMax, locPres,
438     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
439 mlosch 1.6 O bulkMod,
440     I myThid )
441    
442 jmc 1.12 DO j=jMin,jMax
443     DO i=iMin,iMax
444 mlosch 1.6
445     C abbreviations
446 jmc 1.17 t1 = theta(i,j,k,bi,bj)
447 mlosch 1.9 t2 = t1*t1
448     t3 = t2*t1
449 jmc 1.17
450 mlosch 1.9 s1 = salt(i,j,k,bi,bj)
451 jmc 1.14 IF ( s1 .GT. 0. _d 0 ) THEN
452 jmc 1.12 s3o2 = 1.5*SQRT(s1)
453 jmc 1.14 ELSE
454     s1 = 0. _d 0
455     s3o2 = 0. _d 0
456     ENDIF
457 mlosch 1.6
458 jmc 1.12 p1 = locPres(i,j)*SItoBar
459 mlosch 1.6
460     C d(rho)/d(S)
461     C of fresh water at p = 0
462     drhoP0dS = 0. _d 0
463     C of salt water at p = 0
464     drhoP0dS = drhoP0dS
465     & + eosJMDCSw(1)
466 mlosch 1.9 & + eosJMDCSw(2)*t1
467 mlosch 1.6 & + eosJMDCSw(3)*t2
468     & + eosJMDCSw(4)*t3
469 mlosch 1.9 & + eosJMDCSw(5)*t3*t1
470 mlosch 1.6 & + s3o2*(
471     & eosJMDCSw(6)
472 mlosch 1.9 & + eosJMDCSw(7)*t1
473 mlosch 1.6 & + eosJMDCSw(8)*t2
474     & )
475 mlosch 1.9 & + 2*eosJMDCSw(9)*s1
476 mlosch 1.6 C d(bulk modulus)/d(S)
477     C of fresh water at p = 0
478     dKdS = 0. _d 0
479     C of sea water at p = 0
480     dKdSSalt =
481     & eosJMDCKSw(1)
482 mlosch 1.9 & + eosJMDCKSw(2)*t1
483 mlosch 1.6 & + eosJMDCKSw(3)*t2
484     & + eosJMDCKSw(4)*t3
485     & + s3o2*( eosJMDCKSw(5)
486 mlosch 1.9 & + eosJMDCKSw(6)*t1
487 mlosch 1.6 & + eosJMDCKSw(7)*t2
488     & )
489    
490     C of sea water at p
491 jmc 1.17 dKdSPres =
492 mlosch 1.9 & p1*( eosJMDCKP(5)
493     & + eosJMDCKP(6)*t1
494 mlosch 1.6 & + eosJMDCKP(7)*t2
495     & )
496 mlosch 1.9 & + s3o2*p1*eosJMDCKP(8)
497     & + p1*p1*( eosJMDCKP(12)
498     & + eosJMDCKP(13)*t1
499 mlosch 1.6 & + eosJMDCKP(14)*t2
500     & )
501    
502     dKdS = dKdSSalt + dKdSPres
503    
504 jmc 1.17 betaLoc(i,j) =
505 mlosch 1.6 & ( bulkmod(i,j)**2*drhoP0dS
506 mlosch 1.9 & - bulkmod(i,j)*p1*drhoP0dS
507     & - rhoP0(i,j)*p1*dKdS )
508     & /( bulkmod(i,j) - p1 )**2
509 jmc 1.17
510    
511 jmc 1.12 ENDDO
512     ENDDO
513     ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN
514    
515     CALL PRESSURE_FOR_EOS(
516     I bi, bj, iMin, iMax, jMin, jMax, kRef,
517     O locPres,
518 jmc 1.17 I myThid )
519 jmc 1.12
520 jmc 1.17 CALL FIND_RHONUM(
521     I iMin, iMax, jMin, jMax, locPres,
522     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
523     O rhoLoc,
524     I myThid )
525 heimbach 1.15
526 jmc 1.17 CALL FIND_RHODEN(
527     I iMin, iMax, jMin, jMax, locPres,
528     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
529     O rhoDen,
530     I myThid )
531 heimbach 1.13
532 jmc 1.12 DO j=jMin,jMax
533     DO i=iMin,iMax
534 mlosch 1.9 t1 = theta(i,j,k,bi,bj)
535     t2 = t1*t1
536     s1 = salt(i,j,k,bi,bj)
537 jmc 1.14 IF ( s1 .GT. 0. _d 0 ) THEN
538 jmc 1.17 sp5 = SQRT(s1)
539 jmc 1.14 ELSE
540     s1 = 0. _d 0
541     sp5 = 0. _d 0
542     ENDIF
543 jmc 1.17
544 jmc 1.12 p1 = locPres(i,j)*SItodBar
545 mlosch 1.9 p1t1 = p1*t1
546 jmc 1.17
547     dnum_dsalt = eosMDJWFnum(4)
548 mlosch 1.9 & + eosMDJWFnum(5)*t1
549     & + 2.*eosMDJWFnum(6)*s1 + eosMDJWFnum(9)*p1
550 jmc 1.17 dden_dsalt = eosMDJWFden(5)
551     & + t1*( eosMDJWFden(6) + eosMDJWFden(7)*t2 )
552 mlosch 1.9 & + 1.5*sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2)
553    
554 jmc 1.17 betaLoc(i,j) = rhoDen(i,j)*( dnum_dsalt
555 heimbach 1.15 & - (rhoLoc(i,j)*rhoDen(i,j))*dden_dsalt )
556 mlosch 1.9
557 jmc 1.12 ENDDO
558     ENDDO
559 mlosch 1.9
560 mlosch 1.18 ELSEIF ( equationOfState.EQ.'TEOS10' ) THEN
561    
562     CALL PRESSURE_FOR_EOS(
563     I bi, bj, iMin, iMax, jMin, jMax, kRef,
564     O locPres,
565     I myThid )
566    
567     CALL FIND_RHOTEOS(
568     I iMin, iMax, jMin, jMax, locPres,
569     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
570     O rhoLoc, rhoDen,
571     I myThid )
572    
573     DO j=jMin,jMax
574     DO i=iMin,iMax
575     ct = theta(i,j,k,bi,bj)
576     sa = salt(i,j,k,bi,bj)
577     IF ( sa .GT. 0. _d 0 ) THEN
578     sqrtsa = SQRT(sa)
579     ELSE
580     sa = 0. _d 0
581     sqrtsa = 0. _d 0
582     ENDIF
583     p = locPres(i,j)*SItodBar
584    
585     dnum_dsalt = teos(05) + ct*(teos(06) + teos(07)*ct)
586     & + 1.5*sqrtsa*(teos(08)
587     & + ct*(teos(09) + ct*(teos(10) + teos(11)*ct)))
588     & + p*(teos(15) + teos(16)*ct + p*teos(20))
589    
590     dden_dsalt = teos(26)
591     & + ct*(teos(27) + ct*(teos(28) + ct*(teos(29) + teos(30)*ct)))
592     & + 2.*teos(36)*sa
593     & + 1.5*sqrtsa*(teos(31) + ct*(teos(32) + ct*(teos(33)
594     & + ct*(teos(34) + teos(35)*ct))))
595     & + p*(teos(41) + teos(42)*ct + p*teos(46))
596    
597     betaLoc(i,j) = rhoDen(i,j)*( dnum_dsalt
598     & - (rhoLoc(i,j)*rhoDen(i,j))*dden_dsalt )
599    
600     ENDDO
601     ENDDO
602    
603 jmc 1.12 ELSE
604     WRITE(*,*) 'FIND_BETA: equationOfState = ',equationOfState
605     STOP 'FIND_BETA: "equationOfState" has illegal value'
606     ENDIF
607 adcroft 1.1
608 jmc 1.17 RETURN
609 jmc 1.12 END

  ViewVC Help
Powered by ViewVC 1.1.22