My Project
Public Member Functions | List of all members
eqs_of_state::jacket_mcdougall Interface Reference

Public Member Functions

subroutine jacket_mcdougall_1d (MYS, MYT, MYP, MYRHO)
 
subroutine jacket_mcdougall_2d (MYS, MYT, MYP, MYRHO)
 

Detailed Description

Definition at line 104 of file eqs_of_state.f90.

Member Function/Subroutine Documentation

◆ jacket_mcdougall_1d()

subroutine eqs_of_state::jacket_mcdougall::jacket_mcdougall_1d ( real(sp), dimension(:), intent(in)  MYS,
real(sp), dimension(:), intent(in)  MYT,
real(sp), dimension(:), intent(in)  MYP,
real(sp), dimension(:), intent(inout)  MYRHO 
)

Definition at line 571 of file eqs_of_state.f90.

571 
572  !==============================================================================|
573  IMPLICIT NONE
574  REAL(SP), INTENT(IN), DIMENSION(:) :: MYS,MYT,MYP
575  REAL(SP), INTENT(INOUT), DIMENSION(:):: MYRHO
576 
577  REAL(SP) :: TF,SF,sqrtSF,PBAR,TEMP(10),BULK,BULK0,BULK1,BULK2
578  INTEGER :: K,ub1,lb1
579  !==============================================================================|
580  ! Polynomial expansion coefficients for the computation of in situ |
581  ! density via the nonlinear equation of state for seawater as a |
582  ! function of potential temperature, salinity, and pressure (Jackett |
583  ! and McDougall, 1995). |
584  REAL(SP), PARAMETER :: A00 = +1.965933e+04_sp
585  REAL(SP), PARAMETER :: A01 = +1.444304e+02_sp
586  REAL(SP), PARAMETER :: A02 = -1.706103e+00_sp
587  REAL(SP), PARAMETER :: A03 = +9.648704e-03_sp
588  REAL(SP), PARAMETER :: A04 = -4.190253e-05_sp
589  REAL(SP), PARAMETER :: B00 = +5.284855e+01_sp
590  REAL(SP), PARAMETER :: B01 = -3.101089e-01_sp
591  REAL(SP), PARAMETER :: B02 = +6.283263e-03_sp
592  REAL(SP), PARAMETER :: B03 = -5.084188e-05_sp
593  REAL(SP), PARAMETER :: D00 = +3.886640e-01_sp
594  REAL(SP), PARAMETER :: D01 = +9.085835e-03_sp
595  REAL(SP), PARAMETER :: D02 = -4.619924e-04_sp
596  REAL(SP), PARAMETER :: E00 = +3.186519e+00_sp
597  REAL(SP), PARAMETER :: E01 = +2.212276e-02_sp
598  REAL(SP), PARAMETER :: E02 = -2.984642e-04_sp
599  REAL(SP), PARAMETER :: E03 = +1.956415e-06_sp
600  REAL(SP), PARAMETER :: F00 = +6.704388e-03_sp
601  REAL(SP), PARAMETER :: F01 = -1.847318e-04_sp
602  REAL(SP), PARAMETER :: F02 = +2.059331e-07_sp
603  REAL(SP), PARAMETER :: G00 = +1.480266e-04_sp
604  REAL(SP), PARAMETER :: G01 = +2.102898e-04_sp
605  REAL(SP), PARAMETER :: G02 = -1.202016e-05_sp
606  REAL(SP), PARAMETER :: G03 = +1.394680e-07_sp
607  REAL(SP), PARAMETER :: H00 = -2.040237e-06_sp
608  REAL(SP), PARAMETER :: H01 = +6.128773e-08_sp
609  REAL(SP), PARAMETER :: H02 = +6.207323e-10_sp
610 
611  REAL(SP), PARAMETER :: Q00 = +9.99842594e+02_sp
612  REAL(SP), PARAMETER :: Q01 = +6.793952e-02_sp
613  REAL(SP), PARAMETER :: Q02 = -9.095290e-03_sp
614  REAL(SP), PARAMETER :: Q03 = +1.001685e-04_sp
615  REAL(SP), PARAMETER :: Q04 = -1.120083e-06_sp
616  REAL(SP), PARAMETER :: Q05 = +6.536332e-09_sp
617  REAL(SP), PARAMETER :: U00 = +8.24493e-01_sp
618  REAL(SP), PARAMETER :: U01 = -4.08990e-03_sp
619  REAL(SP), PARAMETER :: U02 = +7.64380e-05_sp
620  REAL(SP), PARAMETER :: U03 = -8.24670e-07_sp
621  REAL(SP), PARAMETER :: U04 = +5.38750e-09_sp
622  REAL(SP), PARAMETER :: V00 = -5.72466e-03_sp
623  REAL(SP), PARAMETER :: V01 = +1.02270e-04_sp
624  REAL(SP), PARAMETER :: V02 = -1.65460e-06_sp
625  REAL(SP), PARAMETER :: W00 = +4.8314e-04_sp
626 
627 
628  ! SET DIMS USING MY S
629  lb1 = lbound(mys,1)
630  ub1 = ubound(mys,1)
631 
632  ! CHECK MY T
633  if(lb1 /= lbound(myt,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
634  if(ub1 /= ubound(myt,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
635 
636  ! CHECK MY P
637  if(lb1 /= lbound(myp,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
638  if(ub1 /= ubound(myp,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
639 
640  ! CHECK MY RHO
641  if(lb1 /= lbound(myrho,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
642  if(ub1 /= ubound(myrho,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
643 
644  !
645  ! CALCULATE DENSITY FROM EQUATION OF STATE
646  !
647  DO k=lb1,ub1
648  tf = myt(k)
649  sf = mys(k)
650  sqrtsf = sqrt(sf)
651 
652  pbar = myp(k)
653 
654  ! Compute density (kg/m3) at standard one atmosphere pressure
655  temp(1)=q00+tf*(q01+tf*(q02+tf*(q03+tf*(q04+tf*q05))))
656  temp(2)=u00+tf*(u01+tf*(u02+tf*(u03+tf*u04)))
657  temp(3)=v00+tf*(v01+tf*v02)
658  myrho(k)=temp(1)+sf*(temp(2)+sqrtsf*temp(3)+sf*w00)
659 
660  ! Compute secant bulk modulus (BULK = BULK0 + BULK1*PBAR + BULK2*PBAR*PBAR)
661  temp(4)=a00+tf*(a01+tf*(a02+tf*(a03+tf*a04)))
662  temp(5)=b00+tf*(b01+tf*(b02+tf*b03))
663  temp(6)=d00+tf*(d01+tf*d02)
664  temp(7)=e00+tf*(e01+tf*(e02+tf*e03))
665  temp(8)=f00+tf*(f01+tf*f02)
666  temp(9)=g01+tf*(g02+tf*g03)
667  temp(10)=h00+tf*(h01+tf*h02)
668 
669  bulk0=temp(4)+sf*(temp(5)+sqrtsf*temp(6))
670  bulk1=temp(7)+sf*(temp(8)+sqrtsf*g00)
671  bulk2=temp(9)+sf*temp(10)
672  bulk = bulk0 + pbar * (bulk1 + pbar * bulk2)
673 
674  ! Compute "in situ" density anomaly (kg/m3)
675  myrho(k)=(myrho(k)*bulk)/(bulk-pbar)
676  myrho(k)= myrho(k)-1000.0_sp
677  END DO
678 
679  !
680  ! CALCULATE RHO1
681  !
682  myrho = myrho*1.e-3_sp
683 
684 
subroutine fatal_error(ER1, ER2, ER3, ER4)
Definition: mod_utils.f90:230
Here is the call graph for this function:

◆ jacket_mcdougall_2d()

subroutine eqs_of_state::jacket_mcdougall::jacket_mcdougall_2d ( real(sp), dimension(:,:), intent(in)  MYS,
real(sp), dimension(:,:), intent(in)  MYT,
real(sp), dimension(:,:), intent(in)  MYP,
real(sp), dimension(:,:), intent(inout)  MYRHO 
)

Definition at line 438 of file eqs_of_state.f90.

438 
439  !==============================================================================|
440  IMPLICIT NONE
441  REAL(SP), INTENT(IN), DIMENSION(:,:) :: MYS,MYT,MYP
442  REAL(SP), INTENT(INOUT), DIMENSION(:,:):: MYRHO
443 
444  REAL(SP) :: TF,SF,sqrtSF,PBAR,TEMP(10),BULK,BULK0,BULK1,BULK2
445  INTEGER :: I,K,ub1,lb1,ub2,lb2
446  !==============================================================================|
447  ! Polynomial expansion coefficients for the computation of in situ |
448  ! density via the nonlinear equation of state for seawater as a |
449  ! function of potential temperature, salinity, and pressure (Jackett |
450  ! and McDougall, 1995). |
451  REAL(SP), PARAMETER :: A00 = +1.965933e+04_sp
452  REAL(SP), PARAMETER :: A01 = +1.444304e+02_sp
453  REAL(SP), PARAMETER :: A02 = -1.706103e+00_sp
454  REAL(SP), PARAMETER :: A03 = +9.648704e-03_sp
455  REAL(SP), PARAMETER :: A04 = -4.190253e-05_sp
456  REAL(SP), PARAMETER :: B00 = +5.284855e+01_sp
457  REAL(SP), PARAMETER :: B01 = -3.101089e-01_sp
458  REAL(SP), PARAMETER :: B02 = +6.283263e-03_sp
459  REAL(SP), PARAMETER :: B03 = -5.084188e-05_sp
460  REAL(SP), PARAMETER :: D00 = +3.886640e-01_sp
461  REAL(SP), PARAMETER :: D01 = +9.085835e-03_sp
462  REAL(SP), PARAMETER :: D02 = -4.619924e-04_sp
463  REAL(SP), PARAMETER :: E00 = +3.186519e+00_sp
464  REAL(SP), PARAMETER :: E01 = +2.212276e-02_sp
465  REAL(SP), PARAMETER :: E02 = -2.984642e-04_sp
466  REAL(SP), PARAMETER :: E03 = +1.956415e-06_sp
467  REAL(SP), PARAMETER :: F00 = +6.704388e-03_sp
468  REAL(SP), PARAMETER :: F01 = -1.847318e-04_sp
469  REAL(SP), PARAMETER :: F02 = +2.059331e-07_sp
470  REAL(SP), PARAMETER :: G00 = +1.480266e-04_sp
471  REAL(SP), PARAMETER :: G01 = +2.102898e-04_sp
472  REAL(SP), PARAMETER :: G02 = -1.202016e-05_sp
473  REAL(SP), PARAMETER :: G03 = +1.394680e-07_sp
474  REAL(SP), PARAMETER :: H00 = -2.040237e-06_sp
475  REAL(SP), PARAMETER :: H01 = +6.128773e-08_sp
476  REAL(SP), PARAMETER :: H02 = +6.207323e-10_sp
477 
478  REAL(SP), PARAMETER :: Q00 = +9.99842594e+02_sp
479  REAL(SP), PARAMETER :: Q01 = +6.793952e-02_sp
480  REAL(SP), PARAMETER :: Q02 = -9.095290e-03_sp
481  REAL(SP), PARAMETER :: Q03 = +1.001685e-04_sp
482  REAL(SP), PARAMETER :: Q04 = -1.120083e-06_sp
483  REAL(SP), PARAMETER :: Q05 = +6.536332e-09_sp
484  REAL(SP), PARAMETER :: U00 = +8.24493e-01_sp
485  REAL(SP), PARAMETER :: U01 = -4.08990e-03_sp
486  REAL(SP), PARAMETER :: U02 = +7.64380e-05_sp
487  REAL(SP), PARAMETER :: U03 = -8.24670e-07_sp
488  REAL(SP), PARAMETER :: U04 = +5.38750e-09_sp
489  REAL(SP), PARAMETER :: V00 = -5.72466e-03_sp
490  REAL(SP), PARAMETER :: V01 = +1.02270e-04_sp
491  REAL(SP), PARAMETER :: V02 = -1.65460e-06_sp
492  REAL(SP), PARAMETER :: W00 = +4.8314e-04_sp
493 
494 
495  ! SET DIMS USING MY S
496  lb1 = lbound(mys,1)
497  ub1 = ubound(mys,1)
498 
499  lb2 = lbound(mys,2)
500  ub2 = ubound(mys,2)
501 
502  ! CHECK MY T
503  if(lb1 /= lbound(myt,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
504  if(ub1 /= ubound(myt,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
505 
506  if(lb2 /= lbound(myt,2)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
507  if(ub2 /= ubound(myt,2)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
508 
509  ! CHECK MY P
510  if(lb1 /= lbound(myp,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
511  if(ub1 /= ubound(myp,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
512 
513  if(lb2 /= lbound(myp,2)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
514  if(ub2 /= ubound(myp,2)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
515 
516  ! CHECK MY RHO
517  if(lb1 /= lbound(myrho,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
518  if(ub1 /= ubound(myrho,1)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
519 
520  if(lb2 /= lbound(myrho,2)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
521  if(ub2 /= ubound(myrho,2)) CALL fatal_error("EQS_OF_STATE:: Dimension mismatch!")
522 
523 
524 
525  !
526  ! CALCULATE DENSITY FROM EQUATION OF STATE
527  !
528  DO i=lb1,ub1
529  DO k=lb2,ub2
530  tf = myt(i,k)
531  sf = mys(i,k)
532  sqrtsf = sqrt(sf)
533 
534  pbar = myp(i,k)
535 
536  ! Compute density (kg/m3) at standard one atmosphere pressure
537  temp(1)=q00+tf*(q01+tf*(q02+tf*(q03+tf*(q04+tf*q05))))
538  temp(2)=u00+tf*(u01+tf*(u02+tf*(u03+tf*u04)))
539  temp(3)=v00+tf*(v01+tf*v02)
540  myrho(i,k)=temp(1)+sf*(temp(2)+sqrtsf*temp(3)+sf*w00)
541 
542  ! Compute secant bulk modulus (BULK = BULK0 + BULK1*PBAR + BULK2*PBAR*PBAR)
543  temp(4)=a00+tf*(a01+tf*(a02+tf*(a03+tf*a04)))
544  temp(5)=b00+tf*(b01+tf*(b02+tf*b03))
545  temp(6)=d00+tf*(d01+tf*d02)
546  temp(7)=e00+tf*(e01+tf*(e02+tf*e03))
547  temp(8)=f00+tf*(f01+tf*f02)
548  temp(9)=g01+tf*(g02+tf*g03)
549  temp(10)=h00+tf*(h01+tf*h02)
550 
551  bulk0=temp(4)+sf*(temp(5)+sqrtsf*temp(6))
552  bulk1=temp(7)+sf*(temp(8)+sqrtsf*g00)
553  bulk2=temp(9)+sf*temp(10)
554  bulk = bulk0 + pbar * (bulk1 + pbar * bulk2)
555 
556  ! Compute "in situ" density anomaly (kg/m3)
557  myrho(i,k)=(myrho(i,k)*bulk)/(bulk-pbar)
558  myrho(i,k)= myrho(i,k)-1000.0_sp
559  END DO
560  END DO
561 
562  !
563  ! CALCULATE RHO1
564  !
565  myrho = myrho*1.e-3_sp
566 
567 
subroutine fatal_error(ER1, ER2, ER3, ER4)
Definition: mod_utils.f90:230
Here is the call graph for this function:

The documentation for this interface was generated from the following file: