/[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.24 - (show annotations) (download)
Tue Feb 18 15:25:09 2003 UTC (21 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint50c_post, checkpoint50c_pre, checkpoint48i_post, checkpoint51, checkpoint50, checkpoint50d_post, checkpoint50b_pre, checkpoint51d_post, checkpoint48h_post, checkpoint51b_pre, checkpoint50f_post, checkpoint50a_post, checkpoint50f_pre, checkpoint51b_post, checkpoint51c_post, checkpoint50g_post, checkpoint50h_post, checkpoint50e_pre, checkpoint50i_post, checkpoint50e_post, checkpoint50d_pre, checkpoint51e_post, checkpoint49, checkpoint51f_pre, checkpoint48g_post, checkpoint50b_post, checkpoint51a_post
Changes since 1.23: +176 -162 lines
o compute locally the pressure for use in EOS: UNESCO, JMD95P or MDJWF
o store total Potential in totPhyHyd for diagnostic & EOS funct. of P
o fix restart and overlap Pb when using Z-coord and EOS funct. of P

1 C $Header: /u/gcmpack/MITgcm/model/src/find_rho.F,v 1.23 2002/11/15 03:01:21 heimbach 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 alway calculate the perturbation with respect to the surface level.
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't 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't 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 _RL rfresh, rsalt
244 INTEGER myThid
245
246 C !LOCAL VARIABLES:
247 C == Local variables ==
248 INTEGER i,j
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 s3o2 = s*SQRT(s)
262
263 C density of freshwater at the surface
264 rfresh =
265 & eosJMDCFw(1)
266 & + eosJMDCFw(2)*t
267 & + eosJMDCFw(3)*t2
268 & + eosJMDCFw(4)*t3
269 & + eosJMDCFw(5)*t4
270 & + eosJMDCFw(6)*t4*t
271 C density of sea water at the surface
272 rsalt =
273 & s*(
274 & eosJMDCSw(1)
275 & + eosJMDCSw(2)*t
276 & + eosJMDCSw(3)*t2
277 & + eosJMDCSw(4)*t3
278 & + eosJMDCSw(5)*t4
279 & )
280 & + s3o2*(
281 & eosJMDCSw(6)
282 & + eosJMDCSw(7)*t
283 & + eosJMDCSw(8)*t2
284 & )
285 & + eosJMDCSw(9)*s*s
286
287 rhoP0(i,j) = rfresh + rsalt
288
289 ENDDO
290 ENDDO
291
292 RETURN
293 END
294
295 C !ROUTINE: FIND_BULKMOD
296 C !INTERFACE:
297 SUBROUTINE FIND_BULKMOD(
298 I bi, bj, iMin, iMax, jMin, jMax, k,
299 I locPres, tFld, sFld,
300 O bulkMod,
301 I myThid )
302
303 C !DESCRIPTION: \bv
304 C *==========================================================*
305 C | o SUBROUTINE FIND_BULKMOD
306 C | Calculates the secant bulk modulus K(S,T,p) of a slice
307 C *==========================================================*
308 C |
309 C | k - is the level of Theta/Salt slice
310 C |
311 C *==========================================================*
312 C \ev
313
314 C !USES:
315 IMPLICIT NONE
316 C == Global variables ==
317 #include "SIZE.h"
318 #include "EEPARAMS.h"
319 #include "PARAMS.h"
320 #include "EOS.h"
321
322 C !INPUT/OUTPUT PARAMETERS:
323 C == Routine arguments ==
324 C k :: Level of Theta/Salt slice
325 INTEGER bi,bj,iMin,iMax,jMin,jMax
326 INTEGER k
327 INTEGER kRef
328 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
329 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
330 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
331 _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
332 _RL bMfresh, bMsalt, bMpres
333 INTEGER myThid
334
335 C !LOCAL VARIABLES:
336 C == Local variables ==
337 INTEGER i,j
338 _RL t, t2, t3, t4, s, s3o2, p, p2
339 CEOP
340
341 DO j=jMin,jMax
342 DO i=iMin,iMax
343 C abbreviations
344 t = tFld(i,j,k,bi,bj)
345 t2 = t*t
346 t3 = t2*t
347 t4 = t3*t
348
349 s = sFld(i,j,k,bi,bj)
350 s3o2 = s*SQRT(s)
351 C
352 p = locPres(i,j)*SItoBar
353 p2 = p*p
354 C secant bulk modulus of fresh water at the surface
355 bMfresh =
356 & eosJMDCKFw(1)
357 & + eosJMDCKFw(2)*t
358 & + eosJMDCKFw(3)*t2
359 & + eosJMDCKFw(4)*t3
360 & + eosJMDCKFw(5)*t4
361 C secant bulk modulus of sea water at the surface
362 bMsalt =
363 & s*( eosJMDCKSw(1)
364 & + eosJMDCKSw(2)*t
365 & + eosJMDCKSw(3)*t2
366 & + eosJMDCKSw(4)*t3
367 & )
368 & + s3o2*( eosJMDCKSw(5)
369 & + eosJMDCKSw(6)*t
370 & + eosJMDCKSw(7)*t2
371 & )
372 C secant bulk modulus of sea water at pressure p
373 bMpres =
374 & p*( eosJMDCKP(1)
375 & + eosJMDCKP(2)*t
376 & + eosJMDCKP(3)*t2
377 & + eosJMDCKP(4)*t3
378 & )
379 & + p*s*( eosJMDCKP(5)
380 & + eosJMDCKP(6)*t
381 & + eosJMDCKP(7)*t2
382 & )
383 & + p*s3o2*eosJMDCKP(8)
384 & + p2*( eosJMDCKP(9)
385 & + eosJMDCKP(10)*t
386 & + eosJMDCKP(11)*t2
387 & )
388 & + p2*s*( eosJMDCKP(12)
389 & + eosJMDCKP(13)*t
390 & + eosJMDCKP(14)*t2
391 & )
392
393 bulkMod(i,j) = bMfresh + bMsalt + bMpres
394
395 ENDDO
396 ENDDO
397
398 RETURN
399 END
400
401 CBOP
402 C !ROUTINE: FIND_RHONUM
403 C !INTERFACE:
404 SUBROUTINE FIND_RHONUM(
405 I bi, bj, iMin, iMax, jMin, jMax, k,
406 I locPres, tFld, sFld,
407 O rhoNum,
408 I myThid )
409
410 C !DESCRIPTION: \bv
411 C *==========================================================*
412 C | o SUBROUTINE FIND_RHONUM
413 C | Calculates the numerator of the McDougall et al.
414 C | equation of state
415 C | - the code is more or less a copy of MOM4
416 C *==========================================================*
417 C |
418 C | k - is the level of Theta/Salt slice
419 C |
420 C *==========================================================*
421 C \ev
422
423 C !USES:
424 IMPLICIT NONE
425 C == Global variables ==
426 #include "SIZE.h"
427 #include "EEPARAMS.h"
428 #include "PARAMS.h"
429 #include "EOS.h"
430
431 C !INPUT/OUTPUT PARAMETERS:
432 C == Routine arguments ==
433 C k :: Level of Theta/Salt slice
434 INTEGER bi,bj,iMin,iMax,jMin,jMax
435 INTEGER k
436 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
437 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
438 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
439 _RL rhoNum(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
440 INTEGER myThid
441
442 C !LOCAL VARIABLES:
443 C == Local variables ==
444 INTEGER i,j
445 _RL t1, t2, s1, p1
446 CEOP
447 DO j=jMin,jMax
448 DO i=iMin,iMax
449 C abbreviations
450 t1 = tFld(i,j,k,bi,bj)
451 t2 = t1*t1
452 s1 = sFld(i,j,k,bi,bj)
453
454 p1 = locPres(i,j)*SItodBar
455
456 rhoNum(i,j) = eosMDJWFnum(0)
457 & + t1*(eosMDJWFnum(1)
458 & + t1*(eosMDJWFnum(2) + eosMDJWFnum(3)*t1) )
459 & + s1*(eosMDJWFnum(4)
460 & + eosMDJWFnum(5)*t1 + eosMDJWFnum(6)*s1)
461 & + p1*(eosMDJWFnum(7) + eosMDJWFnum(8)*t2
462 & + eosMDJWFnum(9)*s1
463 & + p1*(eosMDJWFnum(10) + eosMDJWFnum(11)*t2) )
464
465 ENDDO
466 ENDDO
467
468 RETURN
469 end
470
471 CBOP
472 C !ROUTINE: FIND_RHODEN
473 C !INTERFACE:
474 SUBROUTINE FIND_RHODEN(
475 I bi, bj, iMin, iMax, jMin, jMax, k,
476 I locPres, tFld, sFld,
477 O rhoDen,
478 I myThid )
479
480 C !DESCRIPTION: \bv
481 C *==========================================================*
482 C | o SUBROUTINE FIND_RHODEN
483 C | Calculates the denominator of the McDougall et al.
484 C | equation of state
485 C | - the code is more or less a copy of MOM4
486 C *==========================================================*
487 C |
488 C | k - is the level of Theta/Salt slice
489 C |
490 C *==========================================================*
491 C \ev
492
493 C !USES:
494 IMPLICIT NONE
495 C == Global variables ==
496 #include "SIZE.h"
497 #include "EEPARAMS.h"
498 #include "PARAMS.h"
499 #include "EOS.h"
500
501 C !INPUT/OUTPUT PARAMETERS:
502 C == Routine arguments ==
503 C k :: Level of Theta/Salt slice
504 INTEGER bi,bj,iMin,iMax,jMin,jMax
505 INTEGER k
506 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
507 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
508 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
509 _RL rhoDen(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
510 INTEGER myThid
511
512 C !LOCAL VARIABLES:
513 C == Local variables ==
514 INTEGER i,j
515 _RL t1, t2, s1, sp5, p1, p1t1
516 _RL den, epsln
517 parameter ( epsln = 0. _d 0 )
518 CEOP
519 DO j=jMin,jMax
520 DO i=iMin,iMax
521 C abbreviations
522 t1 = tFld(i,j,k,bi,bj)
523 t2 = t1*t1
524 s1 = sFld(i,j,k,bi,bj)
525 sp5 = SQRT(s1)
526
527 p1 = locPres(i,j)*SItodBar
528 p1t1 = p1*t1
529
530 den = eosMDJWFden(0)
531 & + t1*(eosMDJWFden(1)
532 & + t1*(eosMDJWFden(2)
533 & + t1*(eosMDJWFden(3) + t1*eosMDJWFden(4) ) ) )
534 & + s1*(eosMDJWFden(5)
535 & + t1*(eosMDJWFden(6)
536 & + eosMDJWFden(7)*t2)
537 & + sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) )
538 & + p1*(eosMDJWFden(10)
539 & + p1t1*(eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) )
540
541 rhoDen(i,j) = 1.0/(epsln+den)
542
543 ENDDO
544 ENDDO
545
546 RETURN
547 end
548
549 SUBROUTINE find_rho_scalar(
550 I tLoc, sLoc, pLoc,
551 O rhoLoc,
552 I myThid )
553
554 C !DESCRIPTION: \bv
555 C *==========================================================*
556 C | o SUBROUTINE FIND_RHO_SCALAR
557 C | Calculates [rho(S,T,p)-rhoConst]
558 C *==========================================================*
559 C \ev
560
561 C !USES:
562 IMPLICIT NONE
563 C == Global variables ==
564 #include "SIZE.h"
565 #include "EEPARAMS.h"
566 #include "PARAMS.h"
567 #include "EOS.h"
568
569 C !INPUT/OUTPUT PARAMETERS:
570 C == Routine arguments ==
571 _RL sLoc, tLoc, pLoc
572 _RL rhoLoc
573 INTEGER myThid
574
575 C !LOCAL VARIABLES:
576 C == Local variables ==
577
578 _RL t1, t2, t3, t4, s1, s3o2, p1, p2, sp5, p1t1
579 _RL rfresh, rsalt, rhoP0
580 _RL bMfresh, bMsalt, bMpres, BulkMod
581 _RL rhoNum, rhoDen, den, epsln
582 parameter ( epsln = 0. _d 0 )
583
584 character*(max_len_mbuf) msgbuf
585 CEOP
586
587 rhoLoc = 0. _d 0
588 rhoP0 = 0. _d 0
589 bulkMod = 0. _d 0
590 rfresh = 0. _d 0
591 rsalt = 0. _d 0
592 bMfresh = 0. _d 0
593 bMsalt = 0. _d 0
594 bMpres = 0. _d 0
595 rhoNum = 0. _d 0
596 rhoDen = 0. _d 0
597 den = 0. _d 0
598
599 t1 = tLoc
600 t2 = t1*t1
601 t3 = t2*t1
602 t4 = t3*t1
603
604 s1 = sLoc
605 IF ( s1 .LT. 0. _d 0 ) THEN
606 C issue a warning
607 WRITE(*,'(a,i3,a,i3,a,i3,a,e13.5)')
608 & ' FIND_RHO_SCALAR: WARNING, salinity = ', s1
609 s1 = 0. _d 0
610 ENDIF
611
612 IF (equationOfState.EQ.'LINEAR') THEN
613
614 rhoLoc = 0. _d 0
615
616 ELSEIF (equationOfState.EQ.'POLY3') THEN
617
618 C this is not correct, there is a field eosSig0 which should be use here
619 C but I DO not intent to include the reference level in this routine
620 WRITE(*,'(a)')
621 & ' FIND_RHO_SCALAR: for POLY3, the density is not'
622 WRITE(*,'(a)')
623 & ' computed correctly in this routine'
624 rhoLoc = 0. _d 0
625
626 ELSEIF ( equationOfState(1:5).EQ.'JMD95'
627 & .OR. equationOfState.EQ.'UNESCO' ) THEN
628 C nonlinear equation of state in pressure coordinates
629
630 s3o2 = s1*SQRT(s1)
631
632 p1 = pLoc*SItoBar
633 p2 = p1*p1
634
635 C density of freshwater at the surface
636 rfresh =
637 & eosJMDCFw(1)
638 & + eosJMDCFw(2)*t1
639 & + eosJMDCFw(3)*t2
640 & + eosJMDCFw(4)*t3
641 & + eosJMDCFw(5)*t4
642 & + eosJMDCFw(6)*t4*t1
643 C density of sea water at the surface
644 rsalt =
645 & s1*(
646 & eosJMDCSw(1)
647 & + eosJMDCSw(2)*t1
648 & + eosJMDCSw(3)*t2
649 & + eosJMDCSw(4)*t3
650 & + eosJMDCSw(5)*t4
651 & )
652 & + s3o2*(
653 & eosJMDCSw(6)
654 & + eosJMDCSw(7)*t1
655 & + eosJMDCSw(8)*t2
656 & )
657 & + eosJMDCSw(9)*s1*s1
658
659 rhoP0 = rfresh + rsalt
660
661 C secant bulk modulus of fresh water at the surface
662 bMfresh =
663 & eosJMDCKFw(1)
664 & + eosJMDCKFw(2)*t1
665 & + eosJMDCKFw(3)*t2
666 & + eosJMDCKFw(4)*t3
667 & + eosJMDCKFw(5)*t4
668 C secant bulk modulus of sea water at the surface
669 bMsalt =
670 & s1*( eosJMDCKSw(1)
671 & + eosJMDCKSw(2)*t1
672 & + eosJMDCKSw(3)*t2
673 & + eosJMDCKSw(4)*t3
674 & )
675 & + s3o2*( eosJMDCKSw(5)
676 & + eosJMDCKSw(6)*t1
677 & + eosJMDCKSw(7)*t2
678 & )
679 C secant bulk modulus of sea water at pressure p
680 bMpres =
681 & p1*( eosJMDCKP(1)
682 & + eosJMDCKP(2)*t1
683 & + eosJMDCKP(3)*t2
684 & + eosJMDCKP(4)*t3
685 & )
686 & + p1*s1*( eosJMDCKP(5)
687 & + eosJMDCKP(6)*t1
688 & + eosJMDCKP(7)*t2
689 & )
690 & + p1*s3o2*eosJMDCKP(8)
691 & + p2*( eosJMDCKP(9)
692 & + eosJMDCKP(10)*t1
693 & + eosJMDCKP(11)*t2
694 & )
695 & + p2*s1*( eosJMDCKP(12)
696 & + eosJMDCKP(13)*t1
697 & + eosJMDCKP(14)*t2
698 & )
699
700 bulkMod = bMfresh + bMsalt + bMpres
701
702 C density of sea water at pressure p
703 rhoLoc = rhoP0/(1. _d 0 - p1/bulkMod) - rhoConst
704
705 ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN
706
707 sp5 = SQRT(s1)
708
709 p1 = pLoc*SItodBar
710 p1t1 = p1*t1
711
712 rhoNum = eosMDJWFnum(0)
713 & + t1*(eosMDJWFnum(1)
714 & + t1*(eosMDJWFnum(2) + eosMDJWFnum(3)*t1) )
715 & + s1*(eosMDJWFnum(4)
716 & + eosMDJWFnum(5)*t1 + eosMDJWFnum(6)*s1)
717 & + p1*(eosMDJWFnum(7) + eosMDJWFnum(8)*t2
718 & + eosMDJWFnum(9)*s1
719 & + p1*(eosMDJWFnum(10) + eosMDJWFnum(11)*t2) )
720
721
722 den = eosMDJWFden(0)
723 & + t1*(eosMDJWFden(1)
724 & + t1*(eosMDJWFden(2)
725 & + t1*(eosMDJWFden(3) + t1*eosMDJWFden(4) ) ) )
726 & + s1*(eosMDJWFden(5)
727 & + t1*(eosMDJWFden(6)
728 & + eosMDJWFden(7)*t2)
729 & + sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) )
730 & + p1*(eosMDJWFden(10)
731 & + p1t1*(eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) )
732
733 rhoDen = 1.0/(epsln+den)
734
735 rhoLoc = rhoNum*rhoDen - rhoConst
736
737 ELSEIF( equationOfState .EQ. 'IDEALG' ) THEN
738 C
739 ELSE
740 WRITE(msgbuf,'(3A)')
741 & ' FIND_RHO_SCALAR : equationOfState = "',
742 & equationOfState,'"'
743 CALL PRINT_ERROR( msgbuf, mythid )
744 STOP 'ABNORMAL END: S/R FIND_RHO_SCALAR'
745 ENDIF
746
747 RETURN
748 END
749
750 CBOP
751 C !ROUTINE: LOOK_FOR_NEG_SALINITY
752 C !INTERFACE:
753 SUBROUTINE LOOK_FOR_NEG_SALINITY(
754 I bi, bj, iMin, iMax, jMin, jMax, k,
755 I sFld,
756 I myThid )
757
758 C !DESCRIPTION: \bv
759 C *==========================================================*
760 C | o SUBROUTINE LOOK_FOR_NEG_SALINITY
761 C | looks for and fixes negative salinity values
762 C | this is necessary IF the equation of state uses
763 C | the square root of salinity
764 C *==========================================================*
765 C |
766 C | k - is the Salt level
767 C |
768 C *==========================================================*
769 C \ev
770
771 C !USES:
772 IMPLICIT NONE
773 C == Global variables ==
774 #include "SIZE.h"
775 #include "EEPARAMS.h"
776 #include "PARAMS.h"
777 #include "GRID.h"
778
779 C !INPUT/OUTPUT PARAMETERS:
780 C == Routine arguments ==
781 C k :: Level of Theta/Salt slice
782 INTEGER bi,bj,iMin,iMax,jMin,jMax
783 INTEGER k
784 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
785 INTEGER myThid
786
787 C !LOCAL VARIABLES:
788 C == Local variables ==
789 INTEGER i,j, localWarning
790 character*(max_len_mbuf) msgbuf
791 CEOP
792
793 localWarning = 0
794 DO j=jMin,jMax
795 DO i=iMin,iMax
796 C abbreviations
797 IF ( sFld(i,j,k,bi,bj) .LT. 0. _d 0 ) THEN
798 localWarning = localWarning + 1
799 sFld(i,j,k,bi,bj) = 0. _d 0
800 ENDIF
801 ENDDO
802 ENDDO
803 C issue a warning
804 IF ( localWarning .GT. 0 ) THEN
805 WRITE(standardMessageUnit,'(A,A)')
806 & 'S/R LOOK_FOR_NEG_SALINITY: found negative salinity',
807 & 'values and reset them to zero.'
808 WRITE(standardMessageUnit,'(A,I3)')
809 & 'S/R LOOK_FOR_NEG_SALINITY: current level is k = ', k
810 ENDIF
811
812 RETURN
813 END

  ViewVC Help
Powered by ViewVC 1.1.22