| 1 |
C $Header: /u/gcmpack/MITgcm/model/src/find_rho.F,v 1.16 2002/08/07 16:55:52 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, eqn, |
| 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)-Rhonil] 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 | eqn - determines the eqn. of state: |
| 26 |
C | 'LINEAR', 'POLY3', 'UNESCO', 'JMD95Z', 'JMD95P' |
| 27 |
C | |
| 28 |
C *==========================================================* |
| 29 |
C \ev |
| 30 |
|
| 31 |
C !USES: |
| 32 |
implicit none |
| 33 |
C == Global variables == |
| 34 |
#include "SIZE.h" |
| 35 |
#include "EEPARAMS.h" |
| 36 |
#include "PARAMS.h" |
| 37 |
#include "EOS.h" |
| 38 |
#include "GRID.h" |
| 39 |
|
| 40 |
C !INPUT/OUTPUT PARAMETERS: |
| 41 |
C == Routine arguments == |
| 42 |
integer bi,bj,iMin,iMax,jMin,jMax |
| 43 |
integer k ! Level of Theta/Salt slice |
| 44 |
integer kRef ! Pressure reference level |
| 45 |
character*(*) eqn |
| 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 |
| 54 |
_RL refTemp,refSalt,sigRef,tP,sP,deltaSig |
| 55 |
_RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
| 56 |
_RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
| 57 |
character*(max_len_mbuf) msgbuf |
| 58 |
CEOP |
| 59 |
|
| 60 |
#ifdef ALLOW_AUTODIFF_TAMC |
| 61 |
DO j=1-OLy,sNy+OLy |
| 62 |
DO i=1-OLx,sNx+OLx |
| 63 |
rholoc(i,j) = 0. _d 0 |
| 64 |
rhoP0(i,j) = 0. _d 0 |
| 65 |
bulkMod(i,j) = 0. _d 0 |
| 66 |
ENDDO |
| 67 |
ENDDO |
| 68 |
#endif |
| 69 |
|
| 70 |
if (eqn.eq.'LINEAR') then |
| 71 |
|
| 72 |
C ***NOTE*** |
| 73 |
C In the linear EOS, to make the static stability calculation meaningful |
| 74 |
C we alway calculate the perturbation with respect to the surface level. |
| 75 |
C ********** |
| 76 |
refTemp=tRef(kRef) |
| 77 |
refSalt=sRef(kRef) |
| 78 |
|
| 79 |
do j=jMin,jMax |
| 80 |
do i=iMin,iMax |
| 81 |
rholoc(i,j)=rhonil*( |
| 82 |
& sBeta*(sFld(i,j,k,bi,bj)-refSalt) |
| 83 |
& -tAlpha*(tFld(i,j,k,bi,bj)-refTemp) ) |
| 84 |
enddo |
| 85 |
enddo |
| 86 |
|
| 87 |
elseif (eqn.eq.'POLY3') then |
| 88 |
|
| 89 |
refTemp=eosRefT(kRef) |
| 90 |
refSalt=eosRefS(kRef) |
| 91 |
sigRef=eosSig0(kRef) + (1000.-Rhonil) |
| 92 |
|
| 93 |
do j=jMin,jMax |
| 94 |
do i=iMin,iMax |
| 95 |
tP=tFld(i,j,k,bi,bj)-refTemp |
| 96 |
sP=sFld(i,j,k,bi,bj)-refSalt |
| 97 |
#ifdef USE_FACTORIZED_POLY |
| 98 |
deltaSig= |
| 99 |
& (( eosC(9,kRef)*sP + eosC(5,kRef) )*sP + eosC(2,kRef) )*sP |
| 100 |
& + ( ( eosC(6,kRef) |
| 101 |
& *tP |
| 102 |
& +eosC(7,kRef)*sP + eosC(3,kRef) |
| 103 |
& )*tP |
| 104 |
& +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef) |
| 105 |
& )*tP |
| 106 |
#else |
| 107 |
deltaSig= |
| 108 |
& eosC(1,kRef)*tP |
| 109 |
& +eosC(2,kRef) *sP |
| 110 |
& +eosC(3,kRef)*tP*tP |
| 111 |
& +eosC(4,kRef)*tP *sP |
| 112 |
& +eosC(5,kRef) *sP*sP |
| 113 |
& +eosC(6,kRef)*tP*tP*tP |
| 114 |
& +eosC(7,kRef)*tP*tP *sP |
| 115 |
& +eosC(8,kRef)*tP *sP*sP |
| 116 |
& +eosC(9,kRef) *sP*sP*sP |
| 117 |
#endif |
| 118 |
rholoc(i,j)=sigRef+deltaSig |
| 119 |
enddo |
| 120 |
enddo |
| 121 |
|
| 122 |
elseif ( eqn(1:5).eq.'JMD95' .or. eqn.eq.'UNESCO' ) then |
| 123 |
C nonlinear equation of state in pressure coordinates |
| 124 |
|
| 125 |
CALL FIND_RHOP0( |
| 126 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
| 127 |
I tFld, sFld, |
| 128 |
O rhoP0, |
| 129 |
I myThid ) |
| 130 |
|
| 131 |
CALL FIND_BULKMOD( |
| 132 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
| 133 |
I tFld, sFld, |
| 134 |
O bulkMod, |
| 135 |
I myThid ) |
| 136 |
|
| 137 |
do j=jMin,jMax |
| 138 |
do i=iMin,iMax |
| 139 |
|
| 140 |
C density of sea water at pressure p |
| 141 |
rholoc(i,j) = rhoP0(i,j) |
| 142 |
& /(1. _d 0 - |
| 143 |
& pressure(i,j,k,bi,bj)*SItoBar/bulkMod(i,j)) |
| 144 |
& - rhonil |
| 145 |
|
| 146 |
end do |
| 147 |
end do |
| 148 |
|
| 149 |
else |
| 150 |
write(msgbuf,'(3a)') ' FIND_RHO: eqn = "',eqn,'"' |
| 151 |
call print_error( msgbuf, mythid ) |
| 152 |
stop 'ABNORMAL END: S/R FIND_RHO' |
| 153 |
endif |
| 154 |
|
| 155 |
return |
| 156 |
end |
| 157 |
|
| 158 |
CBOP |
| 159 |
C !ROUTINE: FIND_RHOP0 |
| 160 |
C !INTERFACE: |
| 161 |
subroutine FIND_RHOP0( |
| 162 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
| 163 |
I tFld, sFld, |
| 164 |
O rhoP0, |
| 165 |
I myThid ) |
| 166 |
|
| 167 |
C !DESCRIPTION: \bv |
| 168 |
C *==========================================================* |
| 169 |
C | o SUBROUTINE FIND_RHOP0 |
| 170 |
C | Calculates rho(S,T,0) of a slice |
| 171 |
C *==========================================================* |
| 172 |
C | |
| 173 |
C | k - is the surface level |
| 174 |
C | |
| 175 |
C *==========================================================* |
| 176 |
C \ev |
| 177 |
|
| 178 |
C !USES: |
| 179 |
implicit none |
| 180 |
C == Global variables == |
| 181 |
#include "SIZE.h" |
| 182 |
#include "EEPARAMS.h" |
| 183 |
#include "PARAMS.h" |
| 184 |
#include "EOS.h" |
| 185 |
|
| 186 |
C !INPUT/OUTPUT PARAMETERS: |
| 187 |
C == Routine arguments == |
| 188 |
integer bi,bj,iMin,iMax,jMin,jMax |
| 189 |
integer k ! Level of surface slice |
| 190 |
_RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
| 191 |
_RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
| 192 |
_RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
| 193 |
_RL rfresh, rsalt |
| 194 |
integer myThid |
| 195 |
|
| 196 |
C !LOCAL VARIABLES: |
| 197 |
C == Local variables == |
| 198 |
integer i,j |
| 199 |
_RL t, t2, t3, t4, s, s3o2 |
| 200 |
CEOP |
| 201 |
|
| 202 |
do j=jMin,jMax |
| 203 |
do i=iMin,iMax |
| 204 |
C abbreviations |
| 205 |
t = tFld(i,j,k,bi,bj) |
| 206 |
t2 = t*t |
| 207 |
t3 = t2*t |
| 208 |
t4 = t3*t |
| 209 |
|
| 210 |
s = sFld(i,j,k,bi,bj) |
| 211 |
if ( s .lt. 0. _d 0 ) then |
| 212 |
C issue a warning |
| 213 |
write(*,'(a,i3,a,i3,a,i3,a,e13.5)') |
| 214 |
& ' FIND_RHOP0: WARNING, s(',i,',',j,',',k,') = ', s |
| 215 |
s = 0. _d 0 |
| 216 |
end if |
| 217 |
s3o2 = s*sqrt(s) |
| 218 |
|
| 219 |
C density of freshwater at the surface |
| 220 |
rfresh = |
| 221 |
& eosJMDCFw(1) |
| 222 |
& + eosJMDCFw(2)*t |
| 223 |
& + eosJMDCFw(3)*t2 |
| 224 |
& + eosJMDCFw(4)*t3 |
| 225 |
& + eosJMDCFw(5)*t4 |
| 226 |
& + eosJMDCFw(6)*t4*t |
| 227 |
C density of sea water at the surface |
| 228 |
rsalt = |
| 229 |
& s*( |
| 230 |
& eosJMDCSw(1) |
| 231 |
& + eosJMDCSw(2)*t |
| 232 |
& + eosJMDCSw(3)*t2 |
| 233 |
& + eosJMDCSw(4)*t3 |
| 234 |
& + eosJMDCSw(5)*t4 |
| 235 |
& ) |
| 236 |
& + s3o2*( |
| 237 |
& eosJMDCSw(6) |
| 238 |
& + eosJMDCSw(7)*t |
| 239 |
& + eosJMDCSw(8)*t2 |
| 240 |
& ) |
| 241 |
& + eosJMDCSw(9)*s*s |
| 242 |
|
| 243 |
rhoP0(i,j) = rfresh + rsalt |
| 244 |
|
| 245 |
end do |
| 246 |
end do |
| 247 |
|
| 248 |
return |
| 249 |
end |
| 250 |
|
| 251 |
C !ROUTINE: FIND_BULKMOD |
| 252 |
C !INTERFACE: |
| 253 |
subroutine FIND_BULKMOD( |
| 254 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
| 255 |
I tFld, sFld, |
| 256 |
O bulkMod, |
| 257 |
I myThid ) |
| 258 |
|
| 259 |
C !DESCRIPTION: \bv |
| 260 |
C *==========================================================* |
| 261 |
C | o SUBROUTINE FIND_BULKMOD |
| 262 |
C | Calculates the secant bulk modulus K(S,T,p) of a slice |
| 263 |
C *==========================================================* |
| 264 |
C | |
| 265 |
C | k - is the surface level |
| 266 |
C | |
| 267 |
C *==========================================================* |
| 268 |
C \ev |
| 269 |
|
| 270 |
C !USES: |
| 271 |
implicit none |
| 272 |
C == Global variables == |
| 273 |
#include "SIZE.h" |
| 274 |
#include "EEPARAMS.h" |
| 275 |
#include "PARAMS.h" |
| 276 |
#include "EOS.h" |
| 277 |
|
| 278 |
C !INPUT/OUTPUT PARAMETERS: |
| 279 |
C == Routine arguments == |
| 280 |
integer bi,bj,iMin,iMax,jMin,jMax |
| 281 |
integer k ! Level of surface slice |
| 282 |
_RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
| 283 |
_RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
| 284 |
_RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
| 285 |
_RL bMfresh, bMsalt, bMpres |
| 286 |
integer myThid |
| 287 |
|
| 288 |
C !LOCAL VARIABLES: |
| 289 |
C == Local variables == |
| 290 |
integer i,j |
| 291 |
_RL t, t2, t3, t4, s, s3o2, p, p2 |
| 292 |
CEOP |
| 293 |
|
| 294 |
do j=jMin,jMax |
| 295 |
do i=iMin,iMax |
| 296 |
C abbreviations |
| 297 |
t = tFld(i,j,k,bi,bj) |
| 298 |
t2 = t*t |
| 299 |
t3 = t2*t |
| 300 |
t4 = t3*t |
| 301 |
|
| 302 |
s = sFld(i,j,k,bi,bj) |
| 303 |
if ( s .lt. 0. _d 0 ) then |
| 304 |
C issue a warning |
| 305 |
write(*,'(a,i3,a,i3,a,i3,a,e13.5)') |
| 306 |
& ' FIND_BULKMOD: WARNING, s(',i,',',j,',',k,') = ', s |
| 307 |
s = 0. _d 0 |
| 308 |
end if |
| 309 |
s3o2 = s*sqrt(s) |
| 310 |
C |
| 311 |
p = pressure(i,j,k,bi,bj)*SItoBar |
| 312 |
p2 = p*p |
| 313 |
C secant bulk modulus of fresh water at the surface |
| 314 |
bMfresh = |
| 315 |
& eosJMDCKFw(1) |
| 316 |
& + eosJMDCKFw(2)*t |
| 317 |
& + eosJMDCKFw(3)*t2 |
| 318 |
& + eosJMDCKFw(4)*t3 |
| 319 |
& + eosJMDCKFw(5)*t4 |
| 320 |
C secant bulk modulus of sea water at the surface |
| 321 |
bMsalt = |
| 322 |
& s*( eosJMDCKSw(1) |
| 323 |
& + eosJMDCKSw(2)*t |
| 324 |
& + eosJMDCKSw(3)*t2 |
| 325 |
& + eosJMDCKSw(4)*t3 |
| 326 |
& ) |
| 327 |
& + s3o2*( eosJMDCKSw(5) |
| 328 |
& + eosJMDCKSw(6)*t |
| 329 |
& + eosJMDCKSw(7)*t2 |
| 330 |
& ) |
| 331 |
C secant bulk modulus of sea water at pressure p |
| 332 |
bMpres = |
| 333 |
& p*( eosJMDCKP(1) |
| 334 |
& + eosJMDCKP(2)*t |
| 335 |
& + eosJMDCKP(3)*t2 |
| 336 |
& + eosJMDCKP(4)*t3 |
| 337 |
& ) |
| 338 |
& + p*s*( eosJMDCKP(5) |
| 339 |
& + eosJMDCKP(6)*t |
| 340 |
& + eosJMDCKP(7)*t2 |
| 341 |
& ) |
| 342 |
& + p*s3o2*eosJMDCKP(8) |
| 343 |
& + p2*( eosJMDCKP(9) |
| 344 |
& + eosJMDCKP(10)*t |
| 345 |
& + eosJMDCKP(11)*t2 |
| 346 |
& ) |
| 347 |
& + p2*s*( eosJMDCKP(12) |
| 348 |
& + eosJMDCKP(13)*t |
| 349 |
& + eosJMDCKP(14)*t2 |
| 350 |
& ) |
| 351 |
|
| 352 |
bulkMod(i,j) = bMfresh + bMsalt + bMpres |
| 353 |
|
| 354 |
end do |
| 355 |
end do |
| 356 |
|
| 357 |
return |
| 358 |
end |
| 359 |
|