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

Contents of /MITgcm/model/src/find_rho.F

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


Revision 1.33 - (show annotations) (download)
Thu Jan 26 18:57:45 2006 UTC (18 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58b_post, checkpoint59, checkpoint58f_post, checkpoint58d_post, checkpoint58a_post, checkpoint58y_post, checkpoint58t_post, checkpoint58m_post, checkpoint58w_post, checkpoint58o_post, checkpoint58p_post, checkpoint58q_post, checkpoint58e_post, mitgcm_mapl_00, checkpoint58r_post, checkpoint58n_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint58k_post, checkpoint58v_post, checkpoint58l_post, checkpoint58g_post, checkpoint58x_post, checkpoint59j, checkpoint58h_post, checkpoint58j_post, checkpoint58i_post, checkpoint58c_post, checkpoint58u_post, checkpoint58s_post
Changes since 1.32: +25 -15 lines
remove comments about linear EOS and "surface" level.

1 C $Header: /u/gcmpack/MITgcm/model/src/find_rho.F,v 1.32 2005/11/04 01:19:24 jmc Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5 #define USE_FACTORIZED_POLY
6
7 CBOP
8 C !ROUTINE: FIND_RHO
9 C !INTERFACE:
10 SUBROUTINE FIND_RHO(
11 I bi, bj, iMin, iMax, jMin, jMax, k, kRef,
12 I tFld, sFld,
13 O rholoc,
14 I myThid )
15
16 C !DESCRIPTION: \bv
17 C *==========================================================*
18 C | o SUBROUTINE FIND_RHO
19 C | Calculates [rho(S,T,z)-rhoConst] of a slice
20 C *==========================================================*
21 C |
22 C | k - is the Theta/Salt level
23 C | kRef - determines pressure reference level
24 C | (not used in 'LINEAR' mode)
25 C |
26 C *==========================================================*
27 C \ev
28
29 C !USES:
30 IMPLICIT NONE
31 C == Global variables ==
32 #include "SIZE.h"
33 #include "EEPARAMS.h"
34 #include "PARAMS.h"
35 #include "EOS.h"
36 #include "GRID.h"
37
38 C !INPUT/OUTPUT PARAMETERS:
39 C == Routine arguments ==
40 C k :: Level of Theta/Salt slice
41 C kRef :: Pressure reference level
42 INTEGER bi,bj,iMin,iMax,jMin,jMax
43 INTEGER k
44 INTEGER kRef
45 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
46 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
47 _RL rholoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
48 INTEGER myThid
49
50 C !LOCAL VARIABLES:
51 C == Local variables ==
52 INTEGER i,j
53 _RL refTemp,refSalt,sigRef,tP,sP,deltaSig,dRho
54 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
55 _RL rhoP0 (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
56 _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
57 _RL rhoNum (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
58 _RL rhoDen (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
59 CHARACTER*(MAX_LEN_MBUF) msgBuf
60 CEOP
61
62 #ifdef ALLOW_AUTODIFF_TAMC
63 DO j=1-OLy,sNy+OLy
64 DO i=1-OLx,sNx+OLx
65 rholoc(i,j) = 0. _d 0
66 rhoP0(i,j) = 0. _d 0
67 bulkMod(i,j) = 0. _d 0
68 ENDDO
69 ENDDO
70 #endif
71
72 #ifdef CHECK_SALINITY_FOR_NEGATIVE_VALUES
73 CALL LOOK_FOR_NEG_SALINITY( bi, bj, iMin, iMax, jMin, jMax, k,
74 & sFld, myThid )
75 #endif
76
77 IF (equationOfState.EQ.'LINEAR') THEN
78
79 C ***NOTE***
80 C In the linear EOS, to make the static stability calculation meaningful
81 C we use reference temp & salt from level kRef ;
82 C **********
83 refTemp=tRef(kRef)
84 refSalt=sRef(kRef)
85
86 dRho = rhoNil-rhoConst
87
88 DO j=jMin,jMax
89 DO i=iMin,iMax
90 rholoc(i,j)=rhoNil*(
91 & sBeta*(sFld(i,j,k,bi,bj)-refSalt)
92 & -tAlpha*(tFld(i,j,k,bi,bj)-refTemp) )
93 & + dRho
94 ENDDO
95 ENDDO
96
97 ELSEIF (equationOfState.EQ.'POLY3') THEN
98
99 refTemp=eosRefT(kRef)
100 refSalt=eosRefS(kRef)
101 sigRef=eosSig0(kRef) + (1000.-rhoConst)
102
103 DO j=jMin,jMax
104 DO i=iMin,iMax
105 tP=tFld(i,j,k,bi,bj)-refTemp
106 sP=sFld(i,j,k,bi,bj)-refSalt
107 #ifdef USE_FACTORIZED_POLY
108 deltaSig=
109 & (( eosC(9,kRef)*sP + eosC(5,kRef) )*sP + eosC(2,kRef) )*sP
110 & + ( ( eosC(6,kRef)
111 & *tP
112 & +eosC(7,kRef)*sP + eosC(3,kRef)
113 & )*tP
114 & +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef)
115 & )*tP
116 #else
117 deltaSig=
118 & eosC(1,kRef)*tP
119 & +eosC(2,kRef) *sP
120 & +eosC(3,kRef)*tP*tP
121 & +eosC(4,kRef)*tP *sP
122 & +eosC(5,kRef) *sP*sP
123 & +eosC(6,kRef)*tP*tP*tP
124 & +eosC(7,kRef)*tP*tP *sP
125 & +eosC(8,kRef)*tP *sP*sP
126 & +eosC(9,kRef) *sP*sP*sP
127 #endif
128 rholoc(i,j)=sigRef+deltaSig
129 ENDDO
130 ENDDO
131
132 ELSEIF ( equationOfState(1:5).EQ.'JMD95'
133 & .OR. equationOfState.EQ.'UNESCO' ) THEN
134 C nonlinear equation of state in pressure coordinates
135
136 CALL PRESSURE_FOR_EOS(
137 I bi, bj, iMin, iMax, jMin, jMax, kRef,
138 O locPres,
139 I myThid )
140
141 CALL FIND_RHOP0(
142 I bi, bj, iMin, iMax, jMin, jMax, k,
143 I tFld, sFld,
144 O rhoP0,
145 I myThid )
146
147 CALL FIND_BULKMOD(
148 I bi, bj, iMin, iMax, jMin, jMax, k,
149 I locPres, tFld, sFld,
150 O bulkMod,
151 I myThid )
152
153 #ifdef ALLOW_AUTODIFF_TAMC
154 cph can not DO storing here since find_rho is called multiple times;
155 cph additional recomp. should be acceptable
156 cphCADJ STORE rhoP0(:,:) = comlev1_bibj_k , key=kkey , byte=isbyte
157 cphCADJ STORE bulkMod(:,:) = comlev1_bibj_k , key=kkey , byte=isbyte
158 #endif /* ALLOW_AUTODIFF_TAMC */
159 DO j=jMin,jMax
160 DO i=iMin,iMax
161
162 C density of sea water at pressure p
163 rholoc(i,j) = rhoP0(i,j)
164 & /(1. _d 0 -
165 & locPres(i,j)*SItoBar/bulkMod(i,j))
166 & - rhoConst
167
168 ENDDO
169 ENDDO
170
171 ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN
172
173 CALL PRESSURE_FOR_EOS(
174 I bi, bj, iMin, iMax, jMin, jMax, kRef,
175 O locPres,
176 I myThid )
177
178 CALL FIND_RHONUM( bi, bj, iMin, iMax, jMin, jMax, k,
179 & locPres, tFld, sFld, rhoNum, myThid )
180
181 CALL FIND_RHODEN( bi, bj, iMin, iMax, jMin, jMax, k,
182 & locPres, tFld, sFld, rhoDen, myThid )
183 #ifdef ALLOW_AUTODIFF_TAMC
184 cph can not DO storing here since find_rho is called multiple times;
185 cph additional recomp. should be acceptable
186 cphCADJ STORE rhoNum(:,:) = comlev1_bibj_k , key=kkey , byte=isbyte
187 cphCADJ STORE rhoDen(:,:) = comlev1_bibj_k , key=kkey , byte=isbyte
188 #endif /* ALLOW_AUTODIFF_TAMC */
189 DO j=jMin,jMax
190 DO i=iMin,iMax
191 rholoc(i,j) = rhoNum(i,j)*rhoDen(i,j) - rhoConst
192 ENDDO
193 ENDDO
194
195 ELSEIF( equationOfState .EQ. 'IDEALG' ) THEN
196 C
197 ELSE
198 WRITE(msgBuf,'(3a)')
199 & ' FIND_RHO: equationOfState = "',equationOfState,'"'
200 CALL PRINT_ERROR( msgBuf, myThid )
201 STOP 'ABNORMAL END: S/R FIND_RHO'
202 ENDIF
203
204 RETURN
205 END
206
207 CBOP
208 C !ROUTINE: FIND_RHOP0
209 C !INTERFACE:
210 SUBROUTINE FIND_RHOP0(
211 I bi, bj, iMin, iMax, jMin, jMax, k,
212 I tFld, sFld,
213 O rhoP0,
214 I myThid )
215
216 C !DESCRIPTION: \bv
217 C *==========================================================*
218 C | o SUBROUTINE FIND_RHOP0
219 C | Calculates rho(S,T,0) of a slice
220 C *==========================================================*
221 C |
222 C | k - is the surface level
223 C |
224 C *==========================================================*
225 C \ev
226
227 C !USES:
228 IMPLICIT NONE
229 C == Global variables ==
230 #include "SIZE.h"
231 #include "EEPARAMS.h"
232 #include "PARAMS.h"
233 #include "EOS.h"
234
235 C !INPUT/OUTPUT PARAMETERS:
236 C == Routine arguments ==
237 C k :: Level of Theta/Salt slice
238 INTEGER bi,bj,iMin,iMax,jMin,jMax
239 INTEGER k
240 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
241 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
242 _RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
243 INTEGER myThid
244
245 C !LOCAL VARIABLES:
246 C == Local variables ==
247 INTEGER i,j
248 _RL rfresh, rsalt
249 _RL t, t2, t3, t4, s, s3o2
250 CEOP
251
252 DO j=jMin,jMax
253 DO i=iMin,iMax
254 C abbreviations
255 t = tFld(i,j,k,bi,bj)
256 t2 = t*t
257 t3 = t2*t
258 t4 = t3*t
259
260 s = sFld(i,j,k,bi,bj)
261 IF ( s .GT. 0. _d 0 ) THEN
262 s3o2 = s*SQRT(s)
263 ELSE
264 s = 0. _d 0
265 s3o2 = 0. _d 0
266 ENDIF
267
268 C density of freshwater at the surface
269 rfresh =
270 & eosJMDCFw(1)
271 & + eosJMDCFw(2)*t
272 & + eosJMDCFw(3)*t2
273 & + eosJMDCFw(4)*t3
274 & + eosJMDCFw(5)*t4
275 & + eosJMDCFw(6)*t4*t
276 C density of sea water at the surface
277 rsalt =
278 & s*(
279 & eosJMDCSw(1)
280 & + eosJMDCSw(2)*t
281 & + eosJMDCSw(3)*t2
282 & + eosJMDCSw(4)*t3
283 & + eosJMDCSw(5)*t4
284 & )
285 & + s3o2*(
286 & eosJMDCSw(6)
287 & + eosJMDCSw(7)*t
288 & + eosJMDCSw(8)*t2
289 & )
290 & + eosJMDCSw(9)*s*s
291
292 rhoP0(i,j) = rfresh + rsalt
293
294 ENDDO
295 ENDDO
296
297 RETURN
298 END
299
300 C !ROUTINE: FIND_BULKMOD
301 C !INTERFACE:
302 SUBROUTINE FIND_BULKMOD(
303 I bi, bj, iMin, iMax, jMin, jMax, k,
304 I locPres, tFld, sFld,
305 O bulkMod,
306 I myThid )
307
308 C !DESCRIPTION: \bv
309 C *==========================================================*
310 C | o SUBROUTINE FIND_BULKMOD
311 C | Calculates the secant bulk modulus K(S,T,p) of a slice
312 C *==========================================================*
313 C |
314 C | k - is the level of Theta/Salt slice
315 C |
316 C *==========================================================*
317 C \ev
318
319 C !USES:
320 IMPLICIT NONE
321 C == Global variables ==
322 #include "SIZE.h"
323 #include "EEPARAMS.h"
324 #include "PARAMS.h"
325 #include "EOS.h"
326
327 C !INPUT/OUTPUT PARAMETERS:
328 C == Routine arguments ==
329 C k :: Level of Theta/Salt slice
330 INTEGER bi,bj,iMin,iMax,jMin,jMax
331 INTEGER k
332 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
333 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
334 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
335 _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
336 INTEGER myThid
337
338 C !LOCAL VARIABLES:
339 C == Local variables ==
340 INTEGER i,j
341 _RL bMfresh, bMsalt, bMpres
342 _RL t, t2, t3, t4, s, s3o2, p, p2
343 CEOP
344
345 DO j=jMin,jMax
346 DO i=iMin,iMax
347 C abbreviations
348 t = tFld(i,j,k,bi,bj)
349 t2 = t*t
350 t3 = t2*t
351 t4 = t3*t
352
353 s = sFld(i,j,k,bi,bj)
354 IF ( s .GT. 0. _d 0 ) THEN
355 s3o2 = s*SQRT(s)
356 ELSE
357 s = 0. _d 0
358 s3o2 = 0. _d 0
359 ENDIF
360 C
361 p = locPres(i,j)*SItoBar
362 p2 = p*p
363 C secant bulk modulus of fresh water at the surface
364 bMfresh =
365 & eosJMDCKFw(1)
366 & + eosJMDCKFw(2)*t
367 & + eosJMDCKFw(3)*t2
368 & + eosJMDCKFw(4)*t3
369 & + eosJMDCKFw(5)*t4
370 C secant bulk modulus of sea water at the surface
371 bMsalt =
372 & s*( eosJMDCKSw(1)
373 & + eosJMDCKSw(2)*t
374 & + eosJMDCKSw(3)*t2
375 & + eosJMDCKSw(4)*t3
376 & )
377 & + s3o2*( eosJMDCKSw(5)
378 & + eosJMDCKSw(6)*t
379 & + eosJMDCKSw(7)*t2
380 & )
381 C secant bulk modulus of sea water at pressure p
382 bMpres =
383 & p*( eosJMDCKP(1)
384 & + eosJMDCKP(2)*t
385 & + eosJMDCKP(3)*t2
386 & + eosJMDCKP(4)*t3
387 & )
388 & + p*s*( eosJMDCKP(5)
389 & + eosJMDCKP(6)*t
390 & + eosJMDCKP(7)*t2
391 & )
392 & + p*s3o2*eosJMDCKP(8)
393 & + p2*( eosJMDCKP(9)
394 & + eosJMDCKP(10)*t
395 & + eosJMDCKP(11)*t2
396 & )
397 & + p2*s*( eosJMDCKP(12)
398 & + eosJMDCKP(13)*t
399 & + eosJMDCKP(14)*t2
400 & )
401
402 bulkMod(i,j) = bMfresh + bMsalt + bMpres
403
404 ENDDO
405 ENDDO
406
407 RETURN
408 END
409
410 CBOP
411 C !ROUTINE: FIND_RHONUM
412 C !INTERFACE:
413 SUBROUTINE FIND_RHONUM(
414 I bi, bj, iMin, iMax, jMin, jMax, k,
415 I locPres, tFld, sFld,
416 O rhoNum,
417 I myThid )
418
419 C !DESCRIPTION: \bv
420 C *==========================================================*
421 C | o SUBROUTINE FIND_RHONUM
422 C | Calculates the numerator of the McDougall et al.
423 C | equation of state
424 C | - the code is more or less a copy of MOM4
425 C *==========================================================*
426 C |
427 C | k - is the level of Theta/Salt slice
428 C |
429 C *==========================================================*
430 C \ev
431
432 C !USES:
433 IMPLICIT NONE
434 C == Global variables ==
435 #include "SIZE.h"
436 #include "EEPARAMS.h"
437 #include "PARAMS.h"
438 #include "EOS.h"
439
440 C !INPUT/OUTPUT PARAMETERS:
441 C == Routine arguments ==
442 C k :: Level of Theta/Salt slice
443 INTEGER bi,bj,iMin,iMax,jMin,jMax
444 INTEGER k
445 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
446 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
447 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
448 _RL rhoNum(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
449 INTEGER myThid
450
451 C !LOCAL VARIABLES:
452 C == Local variables ==
453 INTEGER i,j
454 _RL t1, t2, s1, p1
455 CEOP
456 DO j=jMin,jMax
457 DO i=iMin,iMax
458 C abbreviations
459 t1 = tFld(i,j,k,bi,bj)
460 t2 = t1*t1
461 s1 = sFld(i,j,k,bi,bj)
462
463 p1 = locPres(i,j)*SItodBar
464
465 rhoNum(i,j) = eosMDJWFnum(0)
466 & + t1*(eosMDJWFnum(1)
467 & + t1*(eosMDJWFnum(2) + eosMDJWFnum(3)*t1) )
468 & + s1*(eosMDJWFnum(4)
469 & + eosMDJWFnum(5)*t1 + eosMDJWFnum(6)*s1)
470 & + p1*(eosMDJWFnum(7) + eosMDJWFnum(8)*t2
471 & + eosMDJWFnum(9)*s1
472 & + p1*(eosMDJWFnum(10) + eosMDJWFnum(11)*t2) )
473
474 ENDDO
475 ENDDO
476
477 RETURN
478 end
479
480 CBOP
481 C !ROUTINE: FIND_RHODEN
482 C !INTERFACE:
483 SUBROUTINE FIND_RHODEN(
484 I bi, bj, iMin, iMax, jMin, jMax, k,
485 I locPres, tFld, sFld,
486 O rhoDen,
487 I myThid )
488
489 C !DESCRIPTION: \bv
490 C *==========================================================*
491 C | o SUBROUTINE FIND_RHODEN
492 C | Calculates the denominator of the McDougall et al.
493 C | equation of state
494 C | - the code is more or less a copy of MOM4
495 C *==========================================================*
496 C |
497 C | k - is the level of Theta/Salt slice
498 C |
499 C *==========================================================*
500 C \ev
501
502 C !USES:
503 IMPLICIT NONE
504 C == Global variables ==
505 #include "SIZE.h"
506 #include "EEPARAMS.h"
507 #include "PARAMS.h"
508 #include "EOS.h"
509
510 C !INPUT/OUTPUT PARAMETERS:
511 C == Routine arguments ==
512 C k :: Level of Theta/Salt slice
513 INTEGER bi,bj,iMin,iMax,jMin,jMax
514 INTEGER k
515 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
516 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
517 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
518 _RL rhoDen(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
519 INTEGER myThid
520
521 C !LOCAL VARIABLES:
522 C == Local variables ==
523 INTEGER i,j
524 _RL t1, t2, s1, sp5, p1, p1t1
525 _RL den, epsln
526 parameter ( epsln = 0. _d 0 )
527 CEOP
528 DO j=jMin,jMax
529 DO i=iMin,iMax
530 C abbreviations
531 t1 = tFld(i,j,k,bi,bj)
532 t2 = t1*t1
533 s1 = sFld(i,j,k,bi,bj)
534 IF ( s1 .GT. 0. _d 0 ) THEN
535 sp5 = SQRT(s1)
536 ELSE
537 s1 = 0. _d 0
538 sp5 = 0. _d 0
539 ENDIF
540
541 p1 = locPres(i,j)*SItodBar
542 p1t1 = p1*t1
543
544 den = eosMDJWFden(0)
545 & + t1*(eosMDJWFden(1)
546 & + t1*(eosMDJWFden(2)
547 & + t1*(eosMDJWFden(3) + t1*eosMDJWFden(4) ) ) )
548 & + s1*(eosMDJWFden(5)
549 & + t1*(eosMDJWFden(6)
550 & + eosMDJWFden(7)*t2)
551 & + sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) )
552 & + p1*(eosMDJWFden(10)
553 & + p1t1*(eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) )
554
555 rhoDen(i,j) = 1.0/(epsln+den)
556
557 ENDDO
558 ENDDO
559
560 RETURN
561 END
562
563 SUBROUTINE FIND_RHO_SCALAR(
564 I tLoc, sLoc, pLoc,
565 O rhoLoc,
566 I myThid )
567
568 C !DESCRIPTION: \bv
569 C *==========================================================*
570 C | o SUBROUTINE FIND_RHO_SCALAR
571 C | Calculates [rho(S,T,p)-rhoConst]
572 C *==========================================================*
573 C \ev
574
575 C !USES:
576 IMPLICIT NONE
577 C == Global variables ==
578 #include "SIZE.h"
579 #include "EEPARAMS.h"
580 #include "PARAMS.h"
581 #include "EOS.h"
582
583 C !INPUT/OUTPUT PARAMETERS:
584 C == Routine arguments ==
585 _RL sLoc, tLoc, pLoc
586 _RL rhoLoc
587 INTEGER myThid
588
589 C !LOCAL VARIABLES:
590 C == Local variables ==
591
592 _RL t1, t2, t3, t4, s1, s3o2, p1, p2, sp5, p1t1
593 _RL rfresh, rsalt, rhoP0
594 _RL bMfresh, bMsalt, bMpres, BulkMod
595 _RL rhoNum, rhoDen, den, epsln
596 parameter ( epsln = 0. _d 0 )
597
598 CHARACTER*(MAX_LEN_MBUF) msgBuf
599 CEOP
600
601 rhoLoc = 0. _d 0
602 rhoP0 = 0. _d 0
603 bulkMod = 0. _d 0
604 rfresh = 0. _d 0
605 rsalt = 0. _d 0
606 bMfresh = 0. _d 0
607 bMsalt = 0. _d 0
608 bMpres = 0. _d 0
609 rhoNum = 0. _d 0
610 rhoDen = 0. _d 0
611 den = 0. _d 0
612
613 t1 = tLoc
614 t2 = t1*t1
615 t3 = t2*t1
616 t4 = t3*t1
617
618 s1 = sLoc
619 IF ( s1 .LT. 0. _d 0 ) THEN
620 C issue a warning
621 WRITE(msgBuf,'(A,E13.5)')
622 & ' FIND_RHO_SCALAR: WARNING, salinity = ', s1
623 CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
624 & SQUEEZE_RIGHT , myThid )
625 s1 = 0. _d 0
626 ENDIF
627
628 IF (equationOfState.EQ.'LINEAR') THEN
629
630 rholoc = rhoNil*(
631 & sBeta *(sLoc-sRef(1))
632 & -tAlpha*(tLoc-tRef(1))
633 & ) + (rhoNil-rhoConst)
634 c rhoLoc = 0. _d 0
635
636 ELSEIF (equationOfState.EQ.'POLY3') THEN
637
638 C this is not correct, there is a field eosSig0 which should be use here
639 C but I DO not intent to include the reference level in this routine
640 WRITE(msgBuf,'(A)')
641 & ' FIND_RHO_SCALAR: for POLY3, the density is not'
642 CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
643 & SQUEEZE_RIGHT , myThid )
644 WRITE(msgBuf,'(A)')
645 & ' computed correctly in this routine'
646 CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
647 & SQUEEZE_RIGHT , myThid )
648 rhoLoc = 0. _d 0
649
650 ELSEIF ( equationOfState(1:5).EQ.'JMD95'
651 & .OR. equationOfState.EQ.'UNESCO' ) THEN
652 C nonlinear equation of state in pressure coordinates
653
654 s3o2 = s1*SQRT(s1)
655
656 p1 = pLoc*SItoBar
657 p2 = p1*p1
658
659 C density of freshwater at the surface
660 rfresh =
661 & eosJMDCFw(1)
662 & + eosJMDCFw(2)*t1
663 & + eosJMDCFw(3)*t2
664 & + eosJMDCFw(4)*t3
665 & + eosJMDCFw(5)*t4
666 & + eosJMDCFw(6)*t4*t1
667 C density of sea water at the surface
668 rsalt =
669 & s1*(
670 & eosJMDCSw(1)
671 & + eosJMDCSw(2)*t1
672 & + eosJMDCSw(3)*t2
673 & + eosJMDCSw(4)*t3
674 & + eosJMDCSw(5)*t4
675 & )
676 & + s3o2*(
677 & eosJMDCSw(6)
678 & + eosJMDCSw(7)*t1
679 & + eosJMDCSw(8)*t2
680 & )
681 & + eosJMDCSw(9)*s1*s1
682
683 rhoP0 = rfresh + rsalt
684
685 C secant bulk modulus of fresh water at the surface
686 bMfresh =
687 & eosJMDCKFw(1)
688 & + eosJMDCKFw(2)*t1
689 & + eosJMDCKFw(3)*t2
690 & + eosJMDCKFw(4)*t3
691 & + eosJMDCKFw(5)*t4
692 C secant bulk modulus of sea water at the surface
693 bMsalt =
694 & s1*( eosJMDCKSw(1)
695 & + eosJMDCKSw(2)*t1
696 & + eosJMDCKSw(3)*t2
697 & + eosJMDCKSw(4)*t3
698 & )
699 & + s3o2*( eosJMDCKSw(5)
700 & + eosJMDCKSw(6)*t1
701 & + eosJMDCKSw(7)*t2
702 & )
703 C secant bulk modulus of sea water at pressure p
704 bMpres =
705 & p1*( eosJMDCKP(1)
706 & + eosJMDCKP(2)*t1
707 & + eosJMDCKP(3)*t2
708 & + eosJMDCKP(4)*t3
709 & )
710 & + p1*s1*( eosJMDCKP(5)
711 & + eosJMDCKP(6)*t1
712 & + eosJMDCKP(7)*t2
713 & )
714 & + p1*s3o2*eosJMDCKP(8)
715 & + p2*( eosJMDCKP(9)
716 & + eosJMDCKP(10)*t1
717 & + eosJMDCKP(11)*t2
718 & )
719 & + p2*s1*( eosJMDCKP(12)
720 & + eosJMDCKP(13)*t1
721 & + eosJMDCKP(14)*t2
722 & )
723
724 bulkMod = bMfresh + bMsalt + bMpres
725
726 C density of sea water at pressure p
727 rhoLoc = rhoP0/(1. _d 0 - p1/bulkMod) - rhoConst
728
729 ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN
730
731 sp5 = SQRT(s1)
732
733 p1 = pLoc*SItodBar
734 p1t1 = p1*t1
735
736 rhoNum = eosMDJWFnum(0)
737 & + t1*(eosMDJWFnum(1)
738 & + t1*(eosMDJWFnum(2) + eosMDJWFnum(3)*t1) )
739 & + s1*(eosMDJWFnum(4)
740 & + eosMDJWFnum(5)*t1 + eosMDJWFnum(6)*s1)
741 & + p1*(eosMDJWFnum(7) + eosMDJWFnum(8)*t2
742 & + eosMDJWFnum(9)*s1
743 & + p1*(eosMDJWFnum(10) + eosMDJWFnum(11)*t2) )
744
745
746 den = eosMDJWFden(0)
747 & + t1*(eosMDJWFden(1)
748 & + t1*(eosMDJWFden(2)
749 & + t1*(eosMDJWFden(3) + t1*eosMDJWFden(4) ) ) )
750 & + s1*(eosMDJWFden(5)
751 & + t1*(eosMDJWFden(6)
752 & + eosMDJWFden(7)*t2)
753 & + sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) )
754 & + p1*(eosMDJWFden(10)
755 & + p1t1*(eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) )
756
757 rhoDen = 1.0/(epsln+den)
758
759 rhoLoc = rhoNum*rhoDen - rhoConst
760
761 ELSEIF( equationOfState .EQ. 'IDEALG' ) THEN
762 C
763 ELSE
764 WRITE(msgBuf,'(3A)')
765 & ' FIND_RHO_SCALAR : equationOfState = "',
766 & equationOfState,'"'
767 CALL PRINT_ERROR( msgBuf, myThid )
768 STOP 'ABNORMAL END: S/R FIND_RHO_SCALAR'
769 ENDIF
770
771 RETURN
772 END
773
774 CBOP
775 C !ROUTINE: LOOK_FOR_NEG_SALINITY
776 C !INTERFACE:
777 SUBROUTINE LOOK_FOR_NEG_SALINITY(
778 I bi, bj, iMin, iMax, jMin, jMax, k,
779 I sFld,
780 I myThid )
781
782 C !DESCRIPTION: \bv
783 C *==========================================================*
784 C | o SUBROUTINE LOOK_FOR_NEG_SALINITY
785 C | looks for and fixes negative salinity values
786 C | this is necessary IF the equation of state uses
787 C | the square root of salinity
788 C *==========================================================*
789 C |
790 C | k - is the Salt level
791 C |
792 C *==========================================================*
793 C \ev
794
795 C !USES:
796 IMPLICIT NONE
797 C == Global variables ==
798 #include "SIZE.h"
799 #include "EEPARAMS.h"
800 #include "PARAMS.h"
801 #include "GRID.h"
802
803 C !INPUT/OUTPUT PARAMETERS:
804 C == Routine arguments ==
805 C k :: Level of Theta/Salt slice
806 INTEGER bi,bj,iMin,iMax,jMin,jMax
807 INTEGER k
808 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
809 INTEGER myThid
810
811 C !LOCAL VARIABLES:
812 C == Local variables ==
813 INTEGER i,j, localWarning
814 c CHARACTER*(MAX_LEN_MBUF) msgBuf
815 CEOP
816
817 localWarning = 0
818 DO j=jMin,jMax
819 DO i=iMin,iMax
820 C abbreviations
821 IF ( sFld(i,j,k,bi,bj) .LT. 0. _d 0 ) THEN
822 localWarning = localWarning + 1
823 sFld(i,j,k,bi,bj) = 0. _d 0
824 ENDIF
825 ENDDO
826 ENDDO
827 C issue a warning
828 IF ( localWarning .GT. 0 ) THEN
829 WRITE(standardMessageUnit,'(A,A)')
830 & 'S/R LOOK_FOR_NEG_SALINITY: found negative salinity',
831 & 'values and reset them to zero.'
832 WRITE(standardMessageUnit,'(A,I3)')
833 & 'S/R LOOK_FOR_NEG_SALINITY: current level is k = ', k
834 ENDIF
835
836 RETURN
837 END

  ViewVC Help
Powered by ViewVC 1.1.22