/[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.22 - (show annotations) (download)
Fri Nov 1 22:00:33 2002 UTC (21 years, 7 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint46n_post, checkpoint46l_post, checkpoint46m_post
Changes since 1.21: +70 -22 lines
 made convective adjustment work with pressure coordinates:
 - changed the direction of k-loop in convective_adjustment.F for the
   case of pressure coordinates (OCEANICP,ATMOSPHERIC buoyancyRelation)
 - adjusted the reference pressure k-index in convective_adjustment.F
 - adjusted the convection condition in convect.F (in analogy to
   calc_ivdc.F)
 - convective_adjustment no longer computes anything on the halos
 - removed the warnings about negative salinity from find_rho.F and
   find_alpha.F; instead the new routine look_for_neg_salinity, called
   at the beginning of find_rho, find_alpha, and find_beta, does a
   check of the entire slice, if CPP-option
   CHECK_SALINITY_FOR_NEGATIVE_VALUES is defined

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

  ViewVC Help
Powered by ViewVC 1.1.22