/[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.36 - (show annotations) (download)
Mon Aug 11 22:33:23 2008 UTC (15 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61f, checkpoint61n, checkpoint61q, checkpoint61e, checkpoint61g, checkpoint61d, checkpoint61c, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61x
Changes since 1.35: +13 -73 lines
- remove S/R "FIND_RHO" (replaced by FIND_RHO_2D).

1 C $Header: /u/gcmpack/MITgcm/model/src/find_rho.F,v 1.35 2008/08/09 01:02:28 jmc Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5 #define USE_FACTORIZED_POLY
6
7 C-- File find_rho.F: Routines to compute density
8 C-- Contents
9 C-- o FIND_RHO_2D
10 C-- o FIND_RHOP0
11 C-- o FIND_BULKMOD
12 C-- o FIND_RHONUM
13 C-- o FIND_RHODEN
14 C-- o LOOK_FOR_NEG_SALINITY
15
16 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
17 CBOP
18 C !ROUTINE: FIND_RHO_2D
19 C !INTERFACE:
20 SUBROUTINE FIND_RHO_2D(
21 I iMin, iMax, jMin, jMax, kRef,
22 I tFld, sFld,
23 O rhoLoc,
24 I k, bi, bj, myThid )
25
26 C !DESCRIPTION: \bv
27 C *==========================================================*
28 C | o SUBROUTINE FIND_RHO_2D
29 C | Calculates [rho(S,T,z)-rhoConst] of a 2-D slice
30 C *==========================================================*
31 C |
32 C | kRef - determines pressure reference level
33 C | (not used in 'LINEAR' mode)
34 C | Note: k is not used ; keep it for debugging.
35 C *==========================================================*
36 C \ev
37
38 C !USES:
39 IMPLICIT NONE
40 C == Global variables ==
41 #include "SIZE.h"
42 #include "EEPARAMS.h"
43 #include "PARAMS.h"
44 #include "EOS.h"
45
46 C !INPUT/OUTPUT PARAMETERS:
47 C == Routine arguments ==
48 C k :: Level of Theta/Salt slice
49 C kRef :: Pressure reference level
50 INTEGER iMin,iMax,jMin,jMax
51 INTEGER kRef
52 _RL tFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
53 _RL sFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
54 _RL rhoLoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
55 INTEGER k, bi, bj
56 INTEGER myThid
57
58 C !LOCAL VARIABLES:
59 C == Local variables ==
60 INTEGER i,j
61 _RL refTemp,refSalt,sigRef,tP,sP,deltaSig,dRho
62 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
63 _RL rhoP0 (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
64 _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
65 _RL rhoNum (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
66 _RL rhoDen (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
67 CHARACTER*(MAX_LEN_MBUF) msgBuf
68 CEOP
69
70 #ifdef ALLOW_AUTODIFF_TAMC
71 DO j=1-OLy,sNy+OLy
72 DO i=1-OLx,sNx+OLx
73 rhoLoc(i,j) = 0. _d 0
74 rhoP0(i,j) = 0. _d 0
75 bulkMod(i,j) = 0. _d 0
76 ENDDO
77 ENDDO
78 #endif
79
80 #ifdef CHECK_SALINITY_FOR_NEGATIVE_VALUES
81 CALL LOOK_FOR_NEG_SALINITY(
82 I iMin, iMax, jMin, jMax,
83 U sFld,
84 I k, bi, bj, myThid )
85 #endif
86
87 IF (equationOfState.EQ.'LINEAR') THEN
88
89 C ***NOTE***
90 C In the linear EOS, to make the static stability calculation meaningful
91 C we use reference temp & salt from level kRef ;
92 C **********
93 refTemp=tRef(kRef)
94 refSalt=sRef(kRef)
95
96 dRho = rhoNil-rhoConst
97
98 DO j=jMin,jMax
99 DO i=iMin,iMax
100 rhoLoc(i,j)=rhoNil*(
101 & sBeta*(sFld(i,j)-refSalt)
102 & -tAlpha*(tFld(i,j)-refTemp) )
103 & + dRho
104 ENDDO
105 ENDDO
106
107 ELSEIF (equationOfState.EQ.'POLY3') THEN
108
109 refTemp=eosRefT(kRef)
110 refSalt=eosRefS(kRef)
111 sigRef=eosSig0(kRef) + (1000.-rhoConst)
112
113 DO j=jMin,jMax
114 DO i=iMin,iMax
115 tP=tFld(i,j)-refTemp
116 sP=sFld(i,j)-refSalt
117 #ifdef USE_FACTORIZED_POLY
118 deltaSig=
119 & (( eosC(9,kRef)*sP + eosC(5,kRef) )*sP + eosC(2,kRef) )*sP
120 & + ( ( eosC(6,kRef)
121 & *tP
122 & +eosC(7,kRef)*sP + eosC(3,kRef)
123 & )*tP
124 & +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef)
125 & )*tP
126 #else
127 deltaSig=
128 & eosC(1,kRef)*tP
129 & +eosC(2,kRef) *sP
130 & +eosC(3,kRef)*tP*tP
131 & +eosC(4,kRef)*tP *sP
132 & +eosC(5,kRef) *sP*sP
133 & +eosC(6,kRef)*tP*tP*tP
134 & +eosC(7,kRef)*tP*tP *sP
135 & +eosC(8,kRef)*tP *sP*sP
136 & +eosC(9,kRef) *sP*sP*sP
137 #endif
138 rhoLoc(i,j)=sigRef+deltaSig
139 ENDDO
140 ENDDO
141
142 ELSEIF ( equationOfState(1:5).EQ.'JMD95'
143 & .OR. equationOfState.EQ.'UNESCO' ) THEN
144 C nonlinear equation of state in pressure coordinates
145
146 CALL PRESSURE_FOR_EOS(
147 I bi, bj, iMin, iMax, jMin, jMax, kRef,
148 O locPres,
149 I myThid )
150
151 CALL FIND_RHOP0(
152 I iMin, iMax, jMin, jMax,
153 I tFld, sFld,
154 O rhoP0,
155 I myThid )
156
157 CALL FIND_BULKMOD(
158 I iMin, iMax, jMin, jMax,
159 I locPres, tFld, sFld,
160 O bulkMod,
161 I myThid )
162
163 c#ifdef ALLOW_AUTODIFF_TAMC
164 cph can not DO storing here since find_rho is called multiple times;
165 cph additional recomp. should be acceptable
166 c#endif /* ALLOW_AUTODIFF_TAMC */
167 DO j=jMin,jMax
168 DO i=iMin,iMax
169
170 C density of sea water at pressure p
171 rhoLoc(i,j) = rhoP0(i,j)
172 & /(1. _d 0 -
173 & locPres(i,j)*SItoBar/bulkMod(i,j) )
174 & - rhoConst
175
176 ENDDO
177 ENDDO
178
179 ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN
180
181 CALL PRESSURE_FOR_EOS(
182 I bi, bj, iMin, iMax, jMin, jMax, kRef,
183 O locPres,
184 I myThid )
185
186 CALL FIND_RHONUM(
187 I iMin, iMax, jMin, jMax,
188 I locPres, tFld, sFld,
189 O rhoNum,
190 I myThid )
191
192 CALL FIND_RHODEN(
193 I iMin, iMax, jMin, jMax,
194 I locPres, tFld, sFld,
195 O rhoDen,
196 I myThid )
197
198 c#ifdef ALLOW_AUTODIFF_TAMC
199 cph can not DO storing here since find_rho is called multiple times;
200 cph additional recomp. should be acceptable
201 c#endif /* ALLOW_AUTODIFF_TAMC */
202 DO j=jMin,jMax
203 DO i=iMin,iMax
204 rhoLoc(i,j) = rhoNum(i,j)*rhoDen(i,j) - rhoConst
205 ENDDO
206 ENDDO
207
208 ELSEIF( equationOfState .EQ. 'IDEALG' ) THEN
209 C
210 ELSE
211 WRITE(msgBuf,'(3a)')
212 & ' FIND_RHO_2D: equationOfState = "',equationOfState,'"'
213 CALL PRINT_ERROR( msgBuf, myThid )
214 STOP 'ABNORMAL END: S/R FIND_RHO_2D'
215 ENDIF
216
217 RETURN
218 END
219
220 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
221 CBOP
222 C !ROUTINE: FIND_RHOP0
223 C !INTERFACE:
224 SUBROUTINE FIND_RHOP0(
225 I iMin, iMax, jMin, jMax,
226 I tFld, sFld,
227 O rhoP0,
228 I myThid )
229
230 C !DESCRIPTION: \bv
231 C *==========================================================*
232 C | o SUBROUTINE FIND_RHOP0
233 C | Calculates rho(S,T,0) of a slice
234 C *==========================================================*
235 C *==========================================================*
236 C \ev
237
238 C !USES:
239 IMPLICIT NONE
240 C == Global variables ==
241 #include "SIZE.h"
242 #include "EEPARAMS.h"
243 #include "PARAMS.h"
244 #include "EOS.h"
245
246 C !INPUT/OUTPUT PARAMETERS:
247 C == Routine arguments ==
248 INTEGER iMin,iMax,jMin,jMax
249 _RL tFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
250 _RL sFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
251 _RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
252 INTEGER myThid
253
254 C !LOCAL VARIABLES:
255 C == Local variables ==
256 INTEGER i,j
257 _RL rfresh, rsalt
258 _RL t, t2, t3, t4, s, s3o2
259 CEOP
260
261 DO j=jMin,jMax
262 DO i=iMin,iMax
263 C abbreviations
264 t = tFld(i,j)
265 t2 = t*t
266 t3 = t2*t
267 t4 = t3*t
268
269 s = sFld(i,j)
270 IF ( s .GT. 0. _d 0 ) THEN
271 s3o2 = s*SQRT(s)
272 ELSE
273 s = 0. _d 0
274 s3o2 = 0. _d 0
275 ENDIF
276
277 C density of freshwater at the surface
278 rfresh =
279 & eosJMDCFw(1)
280 & + eosJMDCFw(2)*t
281 & + eosJMDCFw(3)*t2
282 & + eosJMDCFw(4)*t3
283 & + eosJMDCFw(5)*t4
284 & + eosJMDCFw(6)*t4*t
285 C density of sea water at the surface
286 rsalt =
287 & s*(
288 & eosJMDCSw(1)
289 & + eosJMDCSw(2)*t
290 & + eosJMDCSw(3)*t2
291 & + eosJMDCSw(4)*t3
292 & + eosJMDCSw(5)*t4
293 & )
294 & + s3o2*(
295 & eosJMDCSw(6)
296 & + eosJMDCSw(7)*t
297 & + eosJMDCSw(8)*t2
298 & )
299 & + eosJMDCSw(9)*s*s
300
301 rhoP0(i,j) = rfresh + rsalt
302
303 ENDDO
304 ENDDO
305
306 RETURN
307 END
308
309 C !ROUTINE: FIND_BULKMOD
310 C !INTERFACE:
311 SUBROUTINE FIND_BULKMOD(
312 I iMin, iMax, jMin, jMax,
313 I locPres, tFld, sFld,
314 O bulkMod,
315 I myThid )
316
317 C !DESCRIPTION: \bv
318 C *==========================================================*
319 C | o SUBROUTINE FIND_BULKMOD
320 C | Calculates the secant bulk modulus K(S,T,p) of a slice
321 C *==========================================================*
322 C |
323 C | k - is the level of Theta/Salt slice
324 C |
325 C *==========================================================*
326 C \ev
327
328 C !USES:
329 IMPLICIT NONE
330 C == Global variables ==
331 #include "SIZE.h"
332 #include "EEPARAMS.h"
333 #include "PARAMS.h"
334 #include "EOS.h"
335
336 C !INPUT/OUTPUT PARAMETERS:
337 C == Routine arguments ==
338 INTEGER iMin,iMax,jMin,jMax
339 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
340 _RL tFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
341 _RL sFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
342 _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
343 INTEGER myThid
344
345 C !LOCAL VARIABLES:
346 C == Local variables ==
347 INTEGER i,j
348 _RL bMfresh, bMsalt, bMpres
349 _RL t, t2, t3, t4, s, s3o2, p, p2
350 CEOP
351
352 DO j=jMin,jMax
353 DO i=iMin,iMax
354 C abbreviations
355 t = tFld(i,j)
356 t2 = t*t
357 t3 = t2*t
358 t4 = t3*t
359
360 s = sFld(i,j)
361 IF ( s .GT. 0. _d 0 ) THEN
362 s3o2 = s*SQRT(s)
363 ELSE
364 s = 0. _d 0
365 s3o2 = 0. _d 0
366 ENDIF
367 C
368 p = locPres(i,j)*SItoBar
369 p2 = p*p
370 C secant bulk modulus of fresh water at the surface
371 bMfresh =
372 & eosJMDCKFw(1)
373 & + eosJMDCKFw(2)*t
374 & + eosJMDCKFw(3)*t2
375 & + eosJMDCKFw(4)*t3
376 & + eosJMDCKFw(5)*t4
377 C secant bulk modulus of sea water at the surface
378 bMsalt =
379 & s*( eosJMDCKSw(1)
380 & + eosJMDCKSw(2)*t
381 & + eosJMDCKSw(3)*t2
382 & + eosJMDCKSw(4)*t3
383 & )
384 & + s3o2*( eosJMDCKSw(5)
385 & + eosJMDCKSw(6)*t
386 & + eosJMDCKSw(7)*t2
387 & )
388 C secant bulk modulus of sea water at pressure p
389 bMpres =
390 & p*( eosJMDCKP(1)
391 & + eosJMDCKP(2)*t
392 & + eosJMDCKP(3)*t2
393 & + eosJMDCKP(4)*t3
394 & )
395 & + p*s*( eosJMDCKP(5)
396 & + eosJMDCKP(6)*t
397 & + eosJMDCKP(7)*t2
398 & )
399 & + p*s3o2*eosJMDCKP(8)
400 & + p2*( eosJMDCKP(9)
401 & + eosJMDCKP(10)*t
402 & + eosJMDCKP(11)*t2
403 & )
404 & + p2*s*( eosJMDCKP(12)
405 & + eosJMDCKP(13)*t
406 & + eosJMDCKP(14)*t2
407 & )
408
409 bulkMod(i,j) = bMfresh + bMsalt + bMpres
410
411 ENDDO
412 ENDDO
413
414 RETURN
415 END
416
417 CBOP
418 C !ROUTINE: FIND_RHONUM
419 C !INTERFACE:
420 SUBROUTINE FIND_RHONUM(
421 I iMin, iMax, jMin, jMax,
422 I locPres, tFld, sFld,
423 O rhoNum,
424 I myThid )
425
426 C !DESCRIPTION: \bv
427 C *==========================================================*
428 C | o SUBROUTINE FIND_RHONUM
429 C | Calculates the numerator of the McDougall et al.
430 C | equation of state
431 C | - the code is more or less a copy of MOM4
432 C *==========================================================*
433 C *==========================================================*
434 C \ev
435
436 C !USES:
437 IMPLICIT NONE
438 C == Global variables ==
439 #include "SIZE.h"
440 #include "EEPARAMS.h"
441 #include "PARAMS.h"
442 #include "EOS.h"
443
444 C !INPUT/OUTPUT PARAMETERS:
445 C == Routine arguments ==
446 INTEGER iMin,iMax,jMin,jMax
447 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
448 _RL tFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
449 _RL sFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
450 _RL rhoNum (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
451 INTEGER myThid
452
453 C !LOCAL VARIABLES:
454 C == Local variables ==
455 INTEGER i,j
456 _RL t1, t2, s1, p1
457 CEOP
458 DO j=jMin,jMax
459 DO i=iMin,iMax
460 C abbreviations
461 t1 = tFld(i,j)
462 t2 = t1*t1
463 s1 = sFld(i,j)
464
465 p1 = locPres(i,j)*SItodBar
466
467 rhoNum(i,j) = eosMDJWFnum(0)
468 & + t1*(eosMDJWFnum(1)
469 & + t1*(eosMDJWFnum(2) + eosMDJWFnum(3)*t1) )
470 & + s1*(eosMDJWFnum(4)
471 & + eosMDJWFnum(5)*t1 + eosMDJWFnum(6)*s1)
472 & + p1*(eosMDJWFnum(7) + eosMDJWFnum(8)*t2
473 & + eosMDJWFnum(9)*s1
474 & + p1*(eosMDJWFnum(10) + eosMDJWFnum(11)*t2) )
475
476 ENDDO
477 ENDDO
478
479 RETURN
480 END
481
482 CBOP
483 C !ROUTINE: FIND_RHODEN
484 C !INTERFACE:
485 SUBROUTINE FIND_RHODEN(
486 I iMin, iMax, jMin, jMax,
487 I locPres, tFld, sFld,
488 O rhoDen,
489 I myThid )
490
491 C !DESCRIPTION: \bv
492 C *==========================================================*
493 C | o SUBROUTINE FIND_RHODEN
494 C | Calculates the denominator of the McDougall et al.
495 C | equation of state
496 C | - the code is more or less a copy of MOM4
497 C *==========================================================*
498 C *==========================================================*
499 C \ev
500
501 C !USES:
502 IMPLICIT NONE
503 C == Global variables ==
504 #include "SIZE.h"
505 #include "EEPARAMS.h"
506 #include "PARAMS.h"
507 #include "EOS.h"
508
509 C !INPUT/OUTPUT PARAMETERS:
510 C == Routine arguments ==
511 INTEGER iMin,iMax,jMin,jMax
512 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
513 _RL tFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
514 _RL sFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
515 _RL rhoDen (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
516 INTEGER myThid
517
518 C !LOCAL VARIABLES:
519 C == Local variables ==
520 INTEGER i,j
521 _RL t1, t2, s1, sp5, p1, p1t1
522 _RL den, epsln
523 parameter ( epsln = 0. _d 0 )
524 CEOP
525 DO j=jMin,jMax
526 DO i=iMin,iMax
527 C abbreviations
528 t1 = tFld(i,j)
529 t2 = t1*t1
530 s1 = sFld(i,j)
531 IF ( s1 .GT. 0. _d 0 ) THEN
532 sp5 = SQRT(s1)
533 ELSE
534 s1 = 0. _d 0
535 sp5 = 0. _d 0
536 ENDIF
537
538 p1 = locPres(i,j)*SItodBar
539 p1t1 = p1*t1
540
541 den = eosMDJWFden(0)
542 & + t1*(eosMDJWFden(1)
543 & + t1*(eosMDJWFden(2)
544 & + t1*(eosMDJWFden(3) + t1*eosMDJWFden(4) ) ) )
545 & + s1*(eosMDJWFden(5)
546 & + t1*(eosMDJWFden(6)
547 & + eosMDJWFden(7)*t2)
548 & + sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) )
549 & + p1*(eosMDJWFden(10)
550 & + p1t1*(eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) )
551
552 rhoDen(i,j) = 1.0/(epsln+den)
553
554 ENDDO
555 ENDDO
556
557 RETURN
558 END
559
560 CBOP
561 C !ROUTINE: LOOK_FOR_NEG_SALINITY
562 C !INTERFACE:
563 SUBROUTINE LOOK_FOR_NEG_SALINITY(
564 I iMin, iMax, jMin, jMax,
565 U sFld,
566 I k, bi, bj, myThid )
567
568 C !DESCRIPTION: \bv
569 C *==========================================================*
570 C | o SUBROUTINE LOOK_FOR_NEG_SALINITY
571 C | looks for and fixes negative salinity values
572 C | this is necessary IF the equation of state uses
573 C | the square root of salinity
574 C *==========================================================*
575 C |
576 C | k - is the Salt level
577 C |
578 C *==========================================================*
579 C \ev
580
581 C !USES:
582 IMPLICIT NONE
583 C == Global variables ==
584 #include "SIZE.h"
585 #include "EEPARAMS.h"
586 #include "PARAMS.h"
587
588 C !INPUT/OUTPUT PARAMETERS:
589 C == Routine arguments ==
590 C k :: Level of Salt slice
591 INTEGER bi,bj,iMin,iMax,jMin,jMax,k
592 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
593 INTEGER myThid
594
595 C !LOCAL VARIABLES:
596 C == Local variables ==
597 INTEGER i,j, localWarning
598 c CHARACTER*(MAX_LEN_MBUF) msgBuf
599 CEOP
600
601 localWarning = 0
602 DO j=jMin,jMax
603 DO i=iMin,iMax
604 C abbreviations
605 IF ( sFld(i,j) .LT. 0. _d 0 ) THEN
606 localWarning = localWarning + 1
607 sFld(i,j) = 0. _d 0
608 ENDIF
609 ENDDO
610 ENDDO
611 C issue a warning
612 IF ( localWarning .GT. 0 ) THEN
613 WRITE(standardMessageUnit,'(A,I6,A)')
614 & 'S/R LOOK_FOR_NEG_SALINITY: found',localWarning,
615 & ' negative salinity values and reset them to zero.'
616 WRITE(standardMessageUnit,'(A,I3)')
617 & 'S/R LOOK_FOR_NEG_SALINITY: current level is k = ', k
618 ENDIF
619
620 RETURN
621 END

  ViewVC Help
Powered by ViewVC 1.1.22