/[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.30 - (show annotations) (download)
Mon Oct 10 15:56:34 2005 UTC (18 years, 7 months ago) by baylor
Branch: MAIN
Changes since 1.29: +3 -5 lines
Oops.  Move references out of i,j loop.

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

  ViewVC Help
Powered by ViewVC 1.1.22