6 use mom_eos_linear
, only : calculate_density_linear, calculate_spec_vol_linear
7 use mom_eos_linear
, only : calculate_density_derivs_linear
8 use mom_eos_linear
, only : calculate_specvol_derivs_linear, int_density_dz_linear
9 use mom_eos_linear
, only : calculate_density_second_derivs_linear
10 use mom_eos_linear
, only : calculate_compress_linear, int_spec_vol_dp_linear
11 use mom_eos_wright
, only : calculate_density_wright, calculate_spec_vol_wright
12 use mom_eos_wright
, only : calculate_density_derivs_wright
13 use mom_eos_wright
, only : calculate_specvol_derivs_wright, int_density_dz_wright
14 use mom_eos_wright
, only : calculate_compress_wright, int_spec_vol_dp_wright
15 use mom_eos_wright
, only : calculate_density_second_derivs_wright
16 use mom_eos_unesco
, only : calculate_density_unesco, calculate_spec_vol_unesco
17 use mom_eos_unesco
, only : calculate_density_derivs_unesco, calculate_density_unesco
18 use mom_eos_unesco
, only : calculate_compress_unesco
19 use mom_eos_nemo
, only : calculate_density_nemo
20 use mom_eos_nemo
, only : calculate_density_derivs_nemo, calculate_density_nemo
21 use mom_eos_nemo
, only : calculate_compress_nemo
22 use mom_eos_teos10
, only : calculate_density_teos10, calculate_spec_vol_teos10
23 use mom_eos_teos10
, only : calculate_density_derivs_teos10
24 use mom_eos_teos10
, only : calculate_specvol_derivs_teos10
25 use mom_eos_teos10
, only : calculate_density_second_derivs_teos10
26 use mom_eos_teos10
, only : calculate_compress_teos10
27 use mom_eos_teos10
, only : gsw_sp_from_sr, gsw_pt_from_ct
28 use mom_tfreeze
, only : calculate_tfreeze_linear, calculate_tfreeze_millero
29 use mom_tfreeze
, only : calculate_tfreeze_teos10
30 use mom_error_handler
, only : mom_error, fatal, warning, mom_mesg
31 use mom_file_parser
, only : get_param, log_version, param_file_type
32 use mom_hor_index
, only : hor_index_type
33 use mom_string_functions
, only : uppercase
36 implicit none ;
private 38 #include <MOM_memory.h> 109 integer :: form_of_eos = 0
110 integer :: form_of_tfreeze = 0
114 logical :: compressible = .true.
128 real :: kg_m3_to_r = 1.
129 real :: r_to_kg_m3 = 1.
130 real :: rl2_t2_to_pa = 1.
131 real :: l_t_to_m_s = 1.
166 real,
intent(in) :: T
167 real,
intent(in) :: S
168 real,
intent(in) :: pressure
169 real,
intent(out) :: rho
171 real,
optional,
intent(in) :: rho_ref
172 real,
optional,
intent(in) :: scale
178 if (.not.
associated(eos))
call mom_error(fatal, &
179 "calculate_density_scalar called with an unassociated EOS_type EOS.")
181 p_scale = eos%RL2_T2_to_Pa
183 select case (eos%form_of_EOS)
185 call calculate_density_linear(t, s, p_scale*pressure, rho, &
186 eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS, rho_ref)
188 call calculate_density_unesco(t, s, p_scale*pressure, rho, rho_ref)
190 call calculate_density_wright(t, s, p_scale*pressure, rho, rho_ref)
192 call calculate_density_teos10(t, s, p_scale*pressure, rho, rho_ref)
194 call calculate_density_nemo(t, s, p_scale*pressure, rho, rho_ref)
196 call mom_error(fatal,
"calculate_density_scalar: EOS is not valid.")
199 rho_scale = eos%kg_m3_to_R
200 if (
present(scale)) rho_scale = rho_scale * scale
201 rho = rho_scale * rho
212 real,
intent(in) :: T
213 real,
intent(in) :: S
214 real,
intent(in) :: Tvar
215 real,
intent(in) :: TScov
216 real,
intent(in) :: Svar
217 real,
intent(in) :: pressure
218 real,
intent(out) :: rho
220 real,
optional,
intent(in) :: rho_ref
221 real,
optional,
intent(in) :: scale
224 real :: d2RdTT, d2RdST, d2RdSS, d2RdSp, d2RdTp
228 if (.not.
associated(eos))
call mom_error(fatal, &
229 "calculate_stanley_density_scalar called with an unassociated EOS_type EOS.")
231 p_scale = eos%RL2_T2_to_Pa
233 select case (eos%form_of_EOS)
235 call calculate_density_linear(t, s, p_scale*pressure, rho, &
236 eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS, rho_ref)
237 call calculate_density_second_derivs_linear(t, s, pressure, d2rdss, d2rdst, &
238 d2rdtt, d2rdsp, d2rdtp)
240 call calculate_density_wright(t, s, p_scale*pressure, rho, rho_ref)
241 call calculate_density_second_derivs_wright(t, s, pressure, d2rdss, d2rdst, &
242 d2rdtt, d2rdsp, d2rdtp)
244 call calculate_density_teos10(t, s, p_scale*pressure, rho, rho_ref)
245 call calculate_density_second_derivs_teos10(t, s, pressure, d2rdss, d2rdst, &
246 d2rdtt, d2rdsp, d2rdtp)
248 call mom_error(fatal,
"calculate_stanley_density_scalar: EOS is not valid.")
252 rho = rho + ( 0.5 * d2rdtt * tvar + ( d2rdst * tscov + 0.5 * d2rdss * svar ) )
254 rho_scale = eos%kg_m3_to_R
255 if (
present(scale)) rho_scale = rho_scale * scale
256 rho = rho_scale * rho
263 real,
dimension(:),
intent(in) :: T
264 real,
dimension(:),
intent(in) :: S
265 real,
dimension(:),
intent(in) :: pressure
266 real,
dimension(:),
intent(inout) :: rho
267 integer,
intent(in) :: start
268 integer,
intent(in) :: npts
270 real,
optional,
intent(in) :: rho_ref
271 real,
optional,
intent(in) :: scale
275 if (.not.
associated(eos))
call mom_error(fatal, &
276 "calculate_density_array called with an unassociated EOS_type EOS.")
278 select case (eos%form_of_EOS)
280 call calculate_density_linear(t, s, pressure, rho, start, npts, &
281 eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS, rho_ref)
283 call calculate_density_unesco(t, s, pressure, rho, start, npts, rho_ref)
285 call calculate_density_wright(t, s, pressure, rho, start, npts, rho_ref)
287 call calculate_density_teos10(t, s, pressure, rho, start, npts, rho_ref)
289 call calculate_density_nemo(t, s, pressure, rho, start, npts, rho_ref)
291 call mom_error(fatal,
"calculate_density_array: EOS%form_of_EOS is not valid.")
294 if (
present(scale))
then ;
if (scale /= 1.0)
then ;
do j=start,start+npts-1
295 rho(j) = scale * rho(j)
296 enddo ;
endif ;
endif 305 subroutine calculate_stanley_density_array(T, S, pressure, Tvar, TScov, Svar, rho, start, npts, EOS, rho_ref, scale)
306 real,
dimension(:),
intent(in) :: T
307 real,
dimension(:),
intent(in) :: S
308 real,
dimension(:),
intent(in) :: pressure
309 real,
dimension(:),
intent(in) :: Tvar
310 real,
dimension(:),
intent(in) :: TScov
311 real,
dimension(:),
intent(in) :: Svar
312 real,
dimension(:),
intent(inout) :: rho
313 integer,
intent(in) :: start
314 integer,
intent(in) :: npts
316 real,
optional,
intent(in) :: rho_ref
317 real,
optional,
intent(in) :: scale
320 real,
dimension(size(T)) :: d2RdTT, d2RdST, d2RdSS, d2RdSp, d2RdTp
323 if (.not.
associated(eos))
call mom_error(fatal, &
324 "calculate_density_array called with an unassociated EOS_type EOS.")
326 select case (eos%form_of_EOS)
328 call calculate_density_linear(t, s, pressure, rho, start, npts, &
329 eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS, rho_ref)
330 call calculate_density_second_derivs_linear(t, s, pressure, d2rdss, d2rdst, &
331 d2rdtt, d2rdsp, d2rdtp, start, npts)
333 call calculate_density_wright(t, s, pressure, rho, start, npts, rho_ref)
334 call calculate_density_second_derivs_wright(t, s, pressure, d2rdss, d2rdst, &
335 d2rdtt, d2rdsp, d2rdtp, start, npts)
337 call calculate_density_teos10(t, s, pressure, rho, start, npts, rho_ref)
338 call calculate_density_second_derivs_teos10(t, s, pressure, d2rdss, d2rdst, &
339 d2rdtt, d2rdsp, d2rdtp, start, npts)
341 call mom_error(fatal,
"calculate_stanley_density_array: EOS%form_of_EOS is not valid.")
345 do j=start,start+npts-1
347 + ( 0.5 * d2rdtt(j) * tvar(j) + ( d2rdst(j) * tscov(j) + 0.5 * d2rdss(j) * svar(j) ) )
350 if (
present(scale))
then ;
if (scale /= 1.0)
then ;
do j=start,start+npts-1
351 rho(j) = scale * rho(j)
352 enddo ;
endif ;
endif 360 real,
dimension(:),
intent(in) :: T
361 real,
dimension(:),
intent(in) :: S
362 real,
dimension(:),
intent(in) :: pressure
363 real,
dimension(:),
intent(inout) :: rho
365 integer,
dimension(2),
optional,
intent(in) :: dom
367 real,
optional,
intent(in) :: rho_ref
368 real,
optional,
intent(in) :: scale
374 real :: rho_reference
375 real,
dimension(size(rho)) :: pres
376 integer :: i, is, ie, npts
378 if (.not.
associated(eos))
call mom_error(fatal, &
379 "calculate_density_1d called with an unassociated EOS_type EOS.")
381 if (
present(dom))
then 382 is = dom(1) ; ie = dom(2) ; npts = 1 + ie - is
384 is = 1 ; ie =
size(rho) ; npts = 1 + ie - is
387 p_scale = eos%RL2_T2_to_Pa
388 rho_unscale = eos%R_to_kg_m3
390 if ((p_scale == 1.0) .and. (rho_unscale == 1.0))
then 392 elseif (
present(rho_ref))
then 393 do i=is,ie ; pres(i) = p_scale * pressure(i) ;
enddo 394 rho_reference = rho_unscale*rho_ref
398 do i=is,ie ; pres(i) = p_scale * pressure(i) ;
enddo 402 rho_scale = eos%kg_m3_to_R
403 if (
present(scale)) rho_scale = rho_scale * scale
404 if (rho_scale /= 1.0)
then ;
do i=is,ie
405 rho(i) = rho_scale * rho(i)
416 subroutine calculate_stanley_density_1d(T, S, pressure, Tvar, TScov, Svar, rho, EOS, dom, rho_ref, scale)
417 real,
dimension(:),
intent(in) :: T
418 real,
dimension(:),
intent(in) :: S
419 real,
dimension(:),
intent(in) :: pressure
420 real,
dimension(:),
intent(in) :: Tvar
421 real,
dimension(:),
intent(in) :: TScov
422 real,
dimension(:),
intent(in) :: Svar
423 real,
dimension(:),
intent(inout) :: rho
425 integer,
dimension(2),
optional,
intent(in) :: dom
427 real,
optional,
intent(in) :: rho_ref
428 real,
optional,
intent(in) :: scale
433 real,
dimension(size(rho)) :: pres
434 real,
dimension(size(T)) :: d2RdTT, d2RdST, d2RdSS, d2RdSp, d2RdTp
435 integer :: i, is, ie, npts
437 if (.not.
associated(eos))
call mom_error(fatal, &
438 "calculate_density_1d called with an unassociated EOS_type EOS.")
440 if (
present(dom))
then 441 is = dom(1) ; ie = dom(2) ; npts = 1 + ie - is
443 is = 1 ; ie =
size(rho) ; npts = 1 + ie - is
446 p_scale = eos%RL2_T2_to_Pa
448 pres(i) = p_scale * pressure(i)
451 select case (eos%form_of_EOS)
453 call calculate_density_linear(t, s, pres, rho, 1, npts, &
454 eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS, rho_ref)
455 call calculate_density_second_derivs_linear(t, s, pres, d2rdss, d2rdst, &
456 d2rdtt, d2rdsp, d2rdtp, 1, npts)
458 call calculate_density_wright(t, s, pres, rho, 1, npts, rho_ref)
459 call calculate_density_second_derivs_wright(t, s, pres, d2rdss, d2rdst, &
460 d2rdtt, d2rdsp, d2rdtp, 1, npts)
462 call calculate_density_teos10(t, s, pres, rho, 1, npts, rho_ref)
463 call calculate_density_second_derivs_teos10(t, s, pres, d2rdss, d2rdst, &
464 d2rdtt, d2rdsp, d2rdtp, 1, npts)
466 call mom_error(fatal,
"calculate_stanley_density_scalar: EOS is not valid.")
472 + ( 0.5 * d2rdtt(i) * tvar(i) + ( d2rdst(i) * tscov(i) + 0.5 * d2rdss(i) * svar(i) ) )
475 rho_scale = eos%kg_m3_to_R
476 if (
present(scale)) rho_scale = rho_scale * scale
477 if (rho_scale /= 1.0)
then ;
do i=is,ie
478 rho(i) = rho_scale * rho(i)
486 real,
dimension(:),
intent(in) :: T
487 real,
dimension(:),
intent(in) :: S
488 real,
dimension(:),
intent(in) :: pressure
489 real,
dimension(:),
intent(inout) :: specvol
490 integer,
intent(in) :: start
491 integer,
intent(in) :: npts
493 real,
optional,
intent(in) :: spv_ref
494 real,
optional,
intent(in) :: scale
497 real,
dimension(size(specvol)) :: rho
500 if (.not.
associated(eos))
call mom_error(fatal, &
501 "calculate_spec_vol_array called with an unassociated EOS_type EOS.")
503 select case (eos%form_of_EOS)
505 call calculate_spec_vol_linear(t, s, pressure, specvol, start, npts, &
506 eos%rho_T0_S0, eos%drho_dT, eos%drho_dS, spv_ref)
508 call calculate_spec_vol_unesco(t, s, pressure, specvol, start, npts, spv_ref)
510 call calculate_spec_vol_wright(t, s, pressure, specvol, start, npts, spv_ref)
512 call calculate_spec_vol_teos10(t, s, pressure, specvol, start, npts, spv_ref)
514 call calculate_density_nemo(t, s, pressure, rho, start, npts)
515 if (
present(spv_ref))
then 516 specvol(:) = 1.0 / rho(:) - spv_ref
518 specvol(:) = 1.0 / rho(:)
521 call mom_error(fatal,
"calculate_spec_vol_array: EOS%form_of_EOS is not valid.")
524 if (
present(scale))
then ;
if (scale /= 1.0)
then ;
do j=start,start+npts-1
525 specvol(j) = scale * specvol(j)
526 enddo ;
endif ;
endif 533 real,
intent(in) :: T
534 real,
intent(in) :: S
535 real,
intent(in) :: pressure
536 real,
intent(out) :: specvol
538 real,
optional,
intent(in) :: spv_ref
539 real,
optional,
intent(in) :: scale
542 real,
dimension(1) :: Ta, Sa, pres, spv
543 real :: spv_reference
546 if (.not.
associated(eos))
call mom_error(fatal, &
547 "calc_spec_vol_scalar called with an unassociated EOS_type EOS.")
549 pres(1) = eos%RL2_T2_to_Pa*pressure
550 ta(1) = t ; sa(1) = s
552 if (
present(spv_ref))
then 553 spv_reference = eos%kg_m3_to_R*spv_ref
560 spv_scale = eos%R_to_kg_m3
561 if (
present(scale)) spv_scale = spv_scale * scale
562 if (spv_scale /= 1.0)
then 563 specvol = spv_scale * specvol
570 subroutine calc_spec_vol_1d(T, S, pressure, specvol, EOS, dom, spv_ref, scale)
571 real,
dimension(:),
intent(in) :: T
572 real,
dimension(:),
intent(in) :: S
573 real,
dimension(:),
intent(in) :: pressure
574 real,
dimension(:),
intent(inout) :: specvol
576 integer,
dimension(2),
optional,
intent(in) :: dom
578 real,
optional,
intent(in) :: spv_ref
579 real,
optional,
intent(in) :: scale
583 real,
dimension(size(specvol)) :: pres
587 real :: spv_reference
588 integer :: i, is, ie, npts
590 if (.not.
associated(eos))
call mom_error(fatal, &
591 "calc_spec_vol_1d called with an unassociated EOS_type EOS.")
593 if (
present(dom))
then 594 is = dom(1) ; ie = dom(2) ; npts = 1 + ie - is
596 is = 1 ; ie =
size(specvol) ; npts = 1 + ie - is
599 p_scale = eos%RL2_T2_to_Pa
600 spv_unscale = eos%kg_m3_to_R
602 if ((p_scale == 1.0) .and. (spv_unscale == 1.0))
then 604 elseif (
present(spv_ref))
then 605 do i=is,ie ; pres(i) = p_scale * pressure(i) ;
enddo 606 spv_reference = spv_unscale*spv_ref
610 do i=is,ie ; pres(i) = p_scale * pressure(i) ;
enddo 614 spv_scale = eos%R_to_kg_m3
615 if (
present(scale)) spv_scale = spv_scale * scale
616 if (spv_scale /= 1.0)
then ;
do i=is,ie
617 specvol(i) = spv_scale * specvol(i)
625 real,
intent(in) :: S
626 real,
intent(in) :: pressure
627 real,
intent(out) :: T_fr
630 real,
optional,
intent(in) :: pres_scale
635 if (.not.
associated(eos))
call mom_error(fatal, &
636 "calculate_TFreeze_scalar called with an unassociated EOS_type EOS.")
638 p_scale = 1.0 ;
if (
present(pres_scale)) p_scale = pres_scale
640 select case (eos%form_of_TFreeze)
642 call calculate_tfreeze_linear(s, p_scale*pressure, t_fr, eos%TFr_S0_P0, &
643 eos%dTFr_dS, eos%dTFr_dp)
645 call calculate_tfreeze_millero(s, p_scale*pressure, t_fr)
647 call calculate_tfreeze_teos10(s, p_scale*pressure, t_fr)
649 call mom_error(fatal,
"calculate_TFreeze_scalar: form_of_TFreeze is not valid.")
656 real,
dimension(:),
intent(in) :: S
657 real,
dimension(:),
intent(in) :: pressure
658 real,
dimension(:),
intent(inout) :: T_fr
660 integer,
intent(in) :: start
661 integer,
intent(in) :: npts
663 real,
optional,
intent(in) :: pres_scale
666 real,
dimension(size(pressure)) :: pres
670 if (.not.
associated(eos))
call mom_error(fatal, &
671 "calculate_TFreeze_scalar called with an unassociated EOS_type EOS.")
673 p_scale = 1.0 ;
if (
present(pres_scale)) p_scale = pres_scale
675 if (p_scale == 1.0)
then 676 select case (eos%form_of_TFreeze)
678 call calculate_tfreeze_linear(s, pressure, t_fr, start, npts, &
679 eos%TFr_S0_P0, eos%dTFr_dS, eos%dTFr_dp)
681 call calculate_tfreeze_millero(s, pressure, t_fr, start, npts)
683 call calculate_tfreeze_teos10(s, pressure, t_fr, start, npts)
685 call mom_error(fatal,
"calculate_TFreeze_scalar: form_of_TFreeze is not valid.")
688 do j=start,start+npts-1 ; pres(j) = p_scale * pressure(j) ;
enddo 689 select case (eos%form_of_TFreeze)
691 call calculate_tfreeze_linear(s, pres, t_fr, start, npts, &
692 eos%TFr_S0_P0, eos%dTFr_dS, eos%dTFr_dp)
694 call calculate_tfreeze_millero(s, pres, t_fr, start, npts)
696 call calculate_tfreeze_teos10(s, pres, t_fr, start, npts)
698 call mom_error(fatal,
"calculate_TFreeze_scalar: form_of_TFreeze is not valid.")
706 real,
dimension(:),
intent(in) :: T
707 real,
dimension(:),
intent(in) :: S
708 real,
dimension(:),
intent(in) :: pressure
709 real,
dimension(:),
intent(inout) :: drho_dT
711 real,
dimension(:),
intent(inout) :: drho_dS
713 integer,
intent(in) :: start
714 integer,
intent(in) :: npts
716 real,
optional,
intent(in) :: scale
722 if (.not.
associated(eos))
call mom_error(fatal, &
723 "calculate_density_derivs called with an unassociated EOS_type EOS.")
725 select case (eos%form_of_EOS)
727 call calculate_density_derivs_linear(t, s, pressure, drho_dt, drho_ds, eos%Rho_T0_S0, &
728 eos%dRho_dT, eos%dRho_dS, start, npts)
730 call calculate_density_derivs_unesco(t, s, pressure, drho_dt, drho_ds, start, npts)
732 call calculate_density_derivs_wright(t, s, pressure, drho_dt, drho_ds, start, npts)
734 call calculate_density_derivs_teos10(t, s, pressure, drho_dt, drho_ds, start, npts)
736 call calculate_density_derivs_nemo(t, s, pressure, drho_dt, drho_ds, start, npts)
738 call mom_error(fatal,
"calculate_density_derivs_array: EOS%form_of_EOS is not valid.")
741 if (
present(scale))
then ;
if (scale /= 1.0)
then ;
do j=start,start+npts-1
742 drho_dt(j) = scale * drho_dt(j)
743 drho_ds(j) = scale * drho_ds(j)
744 enddo ;
endif ;
endif 751 real,
dimension(:),
intent(in) :: T
752 real,
dimension(:),
intent(in) :: S
753 real,
dimension(:),
intent(in) :: pressure
754 real,
dimension(:),
intent(inout) :: drho_dT
756 real,
dimension(:),
intent(inout) :: drho_dS
759 integer,
dimension(2),
optional,
intent(in) :: dom
761 real,
optional,
intent(in) :: scale
764 real,
dimension(size(drho_dT)) :: pres
767 integer :: i, is, ie, npts
769 if (.not.
associated(eos))
call mom_error(fatal, &
770 "calculate_density_derivs called with an unassociated EOS_type EOS.")
772 if (
present(dom))
then 773 is = dom(1) ; ie = dom(2) ; npts = 1 + ie - is
775 is = 1 ; ie =
size(drho_dt) ; npts = 1 + ie - is
778 p_scale = eos%RL2_T2_to_Pa
780 if (p_scale == 1.0)
then 783 do i=is,ie ; pres(i) = p_scale * pressure(i) ;
enddo 787 rho_scale = eos%kg_m3_to_R
788 if (
present(scale)) rho_scale = rho_scale * scale
789 if (rho_scale /= 1.0)
then ;
do i=is,ie
790 drho_dt(i) = rho_scale * drho_dt(i)
791 drho_ds(i) = rho_scale * drho_ds(i)
800 real,
intent(in) :: T
801 real,
intent(in) :: S
802 real,
intent(in) :: pressure
803 real,
intent(out) :: drho_dT
805 real,
intent(out) :: drho_dS
808 real,
optional,
intent(in) :: scale
815 if (.not.
associated(eos))
call mom_error(fatal, &
816 "calculate_density_derivs called with an unassociated EOS_type EOS.")
818 p_scale = eos%RL2_T2_to_Pa
820 select case (eos%form_of_EOS)
822 call calculate_density_derivs_linear(t, s, p_scale*pressure, drho_dt, drho_ds, &
823 eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS)
825 call calculate_density_derivs_wright(t, s, p_scale*pressure, drho_dt, drho_ds)
827 call calculate_density_derivs_teos10(t, s, p_scale*pressure, drho_dt, drho_ds)
829 call mom_error(fatal,
"calculate_density_derivs_scalar: EOS%form_of_EOS is not valid.")
832 rho_scale = eos%kg_m3_to_R
833 if (
present(scale)) rho_scale = rho_scale * scale
834 if (rho_scale /= 1.0)
then 835 drho_dt = rho_scale * drho_dt
836 drho_ds = rho_scale * drho_ds
843 drho_dS_dP, drho_dT_dP, start, npts, EOS, scale)
844 real,
dimension(:),
intent(in) :: T
845 real,
dimension(:),
intent(in) :: S
846 real,
dimension(:),
intent(in) :: pressure
847 real,
dimension(:),
intent(inout) :: drho_dS_dS
849 real,
dimension(:),
intent(inout) :: drho_dS_dT
851 real,
dimension(:),
intent(inout) :: drho_dT_dT
853 real,
dimension(:),
intent(inout) :: drho_dS_dP
855 real,
dimension(:),
intent(inout) :: drho_dT_dP
857 integer,
intent(in) :: start
858 integer,
intent(in) :: npts
860 real,
optional,
intent(in) :: scale
863 real,
dimension(size(pressure)) :: pres
869 if (.not.
associated(eos))
call mom_error(fatal, &
870 "calculate_density_derivs called with an unassociated EOS_type EOS.")
872 p_scale = eos%RL2_T2_to_Pa
874 if (p_scale == 1.0)
then 875 select case (eos%form_of_EOS)
877 call calculate_density_second_derivs_linear(t, s, pressure, drho_ds_ds, drho_ds_dt, &
878 drho_dt_dt, drho_ds_dp, drho_dt_dp, start, npts)
880 call calculate_density_second_derivs_wright(t, s, pressure, drho_ds_ds, drho_ds_dt, &
881 drho_dt_dt, drho_ds_dp, drho_dt_dp, start, npts)
883 call calculate_density_second_derivs_teos10(t, s, pressure, drho_ds_ds, drho_ds_dt, &
884 drho_dt_dt, drho_ds_dp, drho_dt_dp, start, npts)
886 call mom_error(fatal,
"calculate_density_derivs: EOS%form_of_EOS is not valid.")
889 do j=start,start+npts-1 ; pres(j) = p_scale * pressure(j) ;
enddo 890 select case (eos%form_of_EOS)
892 call calculate_density_second_derivs_linear(t, s, pres, drho_ds_ds, drho_ds_dt, &
893 drho_dt_dt, drho_ds_dp, drho_dt_dp, start, npts)
895 call calculate_density_second_derivs_wright(t, s, pres, drho_ds_ds, drho_ds_dt, &
896 drho_dt_dt, drho_ds_dp, drho_dt_dp, start, npts)
898 call calculate_density_second_derivs_teos10(t, s, pres, drho_ds_ds, drho_ds_dt, &
899 drho_dt_dt, drho_ds_dp, drho_dt_dp, start, npts)
901 call mom_error(fatal,
"calculate_density_derivs: EOS%form_of_EOS is not valid.")
905 rho_scale = eos%kg_m3_to_R
906 if (
present(scale)) rho_scale = rho_scale * scale
907 if (rho_scale /= 1.0)
then ;
do j=start,start+npts-1
908 drho_ds_ds(j) = rho_scale * drho_ds_ds(j)
909 drho_ds_dt(j) = rho_scale * drho_ds_dt(j)
910 drho_dt_dt(j) = rho_scale * drho_dt_dt(j)
911 drho_ds_dp(j) = rho_scale * drho_ds_dp(j)
912 drho_dt_dp(j) = rho_scale * drho_dt_dp(j)
915 if (p_scale /= 1.0)
then 916 i_p_scale = 1.0 / p_scale
917 do j=start,start+npts-1
918 drho_ds_dp(j) = i_p_scale * drho_ds_dp(j)
919 drho_dt_dp(j) = i_p_scale * drho_dt_dp(j)
927 drho_dS_dP, drho_dT_dP, EOS, scale)
928 real,
intent(in) :: T
929 real,
intent(in) :: S
930 real,
intent(in) :: pressure
931 real,
intent(out) :: drho_dS_dS
933 real,
intent(out) :: drho_dS_dT
935 real,
intent(out) :: drho_dT_dT
937 real,
intent(out) :: drho_dS_dP
939 real,
intent(out) :: drho_dT_dP
942 real,
optional,
intent(in) :: scale
949 if (.not.
associated(eos))
call mom_error(fatal, &
950 "calculate_density_derivs called with an unassociated EOS_type EOS.")
952 p_scale = eos%RL2_T2_to_Pa
954 select case (eos%form_of_EOS)
956 call calculate_density_second_derivs_linear(t, s, p_scale*pressure, drho_ds_ds, drho_ds_dt, &
957 drho_dt_dt, drho_ds_dp, drho_dt_dp)
959 call calculate_density_second_derivs_wright(t, s, p_scale*pressure, drho_ds_ds, drho_ds_dt, &
960 drho_dt_dt, drho_ds_dp, drho_dt_dp)
962 call calculate_density_second_derivs_teos10(t, s, p_scale*pressure, drho_ds_ds, drho_ds_dt, &
963 drho_dt_dt, drho_ds_dp, drho_dt_dp)
965 call mom_error(fatal,
"calculate_density_derivs: EOS%form_of_EOS is not valid.")
968 rho_scale = eos%kg_m3_to_R
969 if (
present(scale)) rho_scale = rho_scale * scale
970 if (rho_scale /= 1.0)
then 971 drho_ds_ds = rho_scale * drho_ds_ds
972 drho_ds_dt = rho_scale * drho_ds_dt
973 drho_dt_dt = rho_scale * drho_dt_dt
974 drho_ds_dp = rho_scale * drho_ds_dp
975 drho_dt_dp = rho_scale * drho_dt_dp
978 if (p_scale /= 1.0)
then 979 i_p_scale = 1.0 / p_scale
980 drho_ds_dp = i_p_scale * drho_ds_dp
981 drho_dt_dp = i_p_scale * drho_dt_dp
988 real,
dimension(:),
intent(in) :: T
989 real,
dimension(:),
intent(in) :: S
990 real,
dimension(:),
intent(in) :: pressure
991 real,
dimension(:),
intent(inout) :: dSV_dT
993 real,
dimension(:),
intent(inout) :: dSV_dS
995 integer,
intent(in) :: start
996 integer,
intent(in) :: npts
1000 real,
dimension(size(T)) :: press
1001 real,
dimension(size(T)) :: rho
1002 real,
dimension(size(T)) :: dRho_dT
1003 real,
dimension(size(T)) :: dRho_dS
1006 if (.not.
associated(eos))
call mom_error(fatal, &
1007 "calculate_spec_vol_derivs_array called with an unassociated EOS_type EOS.")
1009 select case (eos%form_of_EOS)
1011 call calculate_specvol_derivs_linear(t, s, pressure, dsv_dt, dsv_ds, start, &
1012 npts, eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS)
1014 call calculate_density_unesco(t, s, pressure, rho, start, npts)
1015 call calculate_density_derivs_unesco(t, s, pressure, drho_dt, drho_ds, start, npts)
1016 do j=start,start+npts-1
1017 dsv_dt(j) = -drho_dt(j)/(rho(j)**2)
1018 dsv_ds(j) = -drho_ds(j)/(rho(j)**2)
1021 call calculate_specvol_derivs_wright(t, s, pressure, dsv_dt, dsv_ds, start, npts)
1023 call calculate_specvol_derivs_teos10(t, s, pressure, dsv_dt, dsv_ds, start, npts)
1025 call calculate_density_nemo(t, s, pressure, rho, start, npts)
1026 call calculate_density_derivs_nemo(t, s, pressure, drho_dt, drho_ds, start, npts)
1027 do j=start,start+npts-1
1028 dsv_dt(j) = -drho_dt(j)/(rho(j)**2)
1029 dsv_ds(j) = -drho_ds(j)/(rho(j)**2)
1032 call mom_error(fatal,
"calculate_spec_vol_derivs_array: EOS%form_of_EOS is not valid.")
1040 real,
dimension(:),
intent(in) :: T
1041 real,
dimension(:),
intent(in) :: S
1042 real,
dimension(:),
intent(in) :: pressure
1043 real,
dimension(:),
intent(inout) :: dSV_dT
1045 real,
dimension(:),
intent(inout) :: dSV_dS
1048 integer,
dimension(2),
optional,
intent(in) :: dom
1050 real,
optional,
intent(in) :: scale
1054 real,
dimension(size(dSV_dT)) :: press
1057 integer :: i, is, ie, npts
1059 if (.not.
associated(eos))
call mom_error(fatal, &
1060 "calculate_spec_vol_derivs_1d called with an unassociated EOS_type EOS.")
1062 if (
present(dom))
then 1063 is = dom(1) ; ie = dom(2) ; npts = 1 + ie - is
1065 is = 1 ; ie =
size(dsv_dt) ; npts = 1 + ie - is
1067 p_scale = eos%RL2_T2_to_Pa
1069 if (p_scale == 1.0)
then 1072 do i=is,ie ; press(i) = p_scale * pressure(i) ;
enddo 1076 spv_scale = eos%R_to_kg_m3
1077 if (
present(scale)) spv_scale = spv_scale * scale
1078 if (spv_scale /= 1.0)
then ;
do i=is,ie
1079 dsv_dt(i) = spv_scale * dsv_dt(i)
1080 dsv_ds(i) = spv_scale * dsv_ds(i)
1089 real,
dimension(:),
intent(in) :: T
1090 real,
dimension(:),
intent(in) :: S
1091 real,
dimension(:),
intent(in) :: press
1092 real,
dimension(:),
intent(inout) :: rho
1093 real,
dimension(:),
intent(inout) :: drho_dp
1096 integer,
intent(in) :: start
1097 integer,
intent(in) :: npts
1101 real,
dimension(size(press)) :: pressure
1102 integer :: i, is, ie
1104 if (.not.
associated(eos))
call mom_error(fatal, &
1105 "calculate_compress called with an unassociated EOS_type EOS.")
1107 is = start ; ie = is + npts - 1
1108 do i=is,ie ; pressure(i) = eos%RL2_T2_to_Pa * press(i) ;
enddo 1110 select case (eos%form_of_EOS)
1112 call calculate_compress_linear(t, s, pressure, rho, drho_dp, start, npts, &
1113 eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS)
1115 call calculate_compress_unesco(t, s, pressure, rho, drho_dp, start, npts)
1117 call calculate_compress_wright(t, s, pressure, rho, drho_dp, start, npts)
1119 call calculate_compress_teos10(t, s, pressure, rho, drho_dp, start, npts)
1121 call calculate_compress_nemo(t, s, pressure, rho, drho_dp, start, npts)
1123 call mom_error(fatal,
"calculate_compress: EOS%form_of_EOS is not valid.")
1126 if (eos%kg_m3_to_R /= 1.0)
then ;
do i=is,ie
1127 rho(i) = eos%kg_m3_to_R * rho(i)
1129 if (eos%L_T_to_m_s /= 1.0)
then ;
do i=is,ie
1130 drho_dp(i) = eos%L_T_to_m_s**2 * drho_dp(i)
1139 real,
intent(in) :: T
1140 real,
intent(in) :: S
1141 real,
intent(in) :: pressure
1142 real,
intent(out) :: rho
1143 real,
intent(out) :: drho_dp
1148 real,
dimension(1) :: Ta, Sa, pa, rhoa, drho_dpa
1150 if (.not.
associated(eos))
call mom_error(fatal, &
1151 "calculate_compress called with an unassociated EOS_type EOS.")
1152 ta(1) = t ; sa(1) = s; pa(1) = pressure
1155 rho = rhoa(1) ; drho_dp = drho_dpa(1)
1163 type(hor_index_type),
intent(in) :: HI
1164 integer,
optional,
intent(in) :: halo
1165 integer,
dimension(2) :: EOSdom
1171 halo_sz = 0 ;
if (
present(halo)) halo_sz = halo
1173 eosdom(1) = hi%isc - (hi%isd-1) - halo_sz
1174 eosdom(2) = hi%iec - (hi%isd-1) + halo_sz
1186 dza, intp_dza, intx_dza, inty_dza, halo_size, &
1187 bathyP, dP_tiny, useMassWghtInterp)
1188 type(hor_index_type),
intent(in) :: HI
1189 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1191 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1193 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1195 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1197 real,
intent(in) :: alpha_ref
1202 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1203 intent(inout) :: dza
1205 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1206 optional,
intent(inout) :: intp_dza
1209 real,
dimension(HI%IsdB:HI%IedB,HI%jsd:HI%jed), &
1210 optional,
intent(inout) :: intx_dza
1213 real,
dimension(HI%isd:HI%ied,HI%JsdB:HI%JedB), &
1214 optional,
intent(inout) :: inty_dza
1217 integer,
optional,
intent(in) :: halo_size
1218 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1219 optional,
intent(in) :: bathyP
1220 real,
optional,
intent(in) :: dP_tiny
1222 logical,
optional,
intent(in) :: useMassWghtInterp
1229 if (.not.
associated(eos))
call mom_error(fatal, &
1230 "int_specific_vol_dp called with an unassociated EOS_type EOS.")
1234 if (eos%EOS_quadrature)
call mom_error(fatal,
"EOS_quadrature is set!")
1236 select case (eos%form_of_EOS)
1238 call int_spec_vol_dp_linear(t, s, p_t, p_b, alpha_ref, hi, eos%kg_m3_to_R*eos%Rho_T0_S0, &
1239 eos%kg_m3_to_R*eos%dRho_dT, eos%kg_m3_to_R*eos%dRho_dS, dza, &
1240 intp_dza, intx_dza, inty_dza, halo_size, &
1241 bathyp, dp_tiny, usemasswghtinterp)
1243 call int_spec_vol_dp_wright(t, s, p_t, p_b, alpha_ref, hi, dza, intp_dza, intx_dza, &
1244 inty_dza, halo_size, bathyp, dp_tiny, usemasswghtinterp, &
1245 sv_scale=eos%R_to_kg_m3, pres_scale=eos%RL2_T2_to_Pa)
1247 call mom_error(fatal,
"No analytic integration option is available with this EOS!")
1255 subroutine analytic_int_density_dz(T, S, z_t, z_b, rho_ref, rho_0, G_e, HI, EOS, dpa, &
1256 intz_dpa, intx_dpa, inty_dpa, bathyT, dz_neglect, useMassWghtInterp)
1257 type(hor_index_type),
intent(in) :: HI
1258 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1260 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1262 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1264 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1266 real,
intent(in) :: rho_ref
1269 real,
intent(in) :: rho_0
1272 real,
intent(in) :: G_e
1275 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1276 intent(inout) :: dpa
1278 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1279 optional,
intent(inout) :: intz_dpa
1282 real,
dimension(HI%IsdB:HI%IedB,HI%jsd:HI%jed), &
1283 optional,
intent(inout) :: intx_dpa
1286 real,
dimension(HI%isd:HI%ied,HI%JsdB:HI%JedB), &
1287 optional,
intent(inout) :: inty_dpa
1290 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed), &
1291 optional,
intent(in) :: bathyT
1292 real,
optional,
intent(in) :: dz_neglect
1293 logical,
optional,
intent(in) :: useMassWghtInterp
1300 if (.not.
associated(eos))
call mom_error(fatal, &
1301 "int_density_dz called with an unassociated EOS_type EOS.")
1305 if (eos%EOS_quadrature)
call mom_error(fatal,
"EOS_quadrature is set!")
1307 select case (eos%form_of_EOS)
1309 rho_scale = eos%kg_m3_to_R
1310 if (rho_scale /= 1.0)
then 1311 call int_density_dz_linear(t, s, z_t, z_b, rho_ref, rho_0, g_e, hi, &
1312 rho_scale*eos%Rho_T0_S0, rho_scale*eos%dRho_dT, rho_scale*eos%dRho_dS, &
1313 dpa, intz_dpa, intx_dpa, inty_dpa, bathyt, dz_neglect, usemasswghtinterp)
1315 call int_density_dz_linear(t, s, z_t, z_b, rho_ref, rho_0, g_e, hi, &
1316 eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS, &
1317 dpa, intz_dpa, intx_dpa, inty_dpa, bathyt, dz_neglect, usemasswghtinterp)
1320 rho_scale = eos%kg_m3_to_R
1321 pres_scale = eos%RL2_T2_to_Pa
1322 if ((rho_scale /= 1.0) .or. (pres_scale /= 1.0))
then 1323 call int_density_dz_wright(t, s, z_t, z_b, rho_ref, rho_0, g_e, hi, &
1324 dpa, intz_dpa, intx_dpa, inty_dpa, bathyt, &
1325 dz_neglect, usemasswghtinterp, rho_scale, pres_scale)
1327 call int_density_dz_wright(t, s, z_t, z_b, rho_ref, rho_0, g_e, hi, &
1328 dpa, intz_dpa, intx_dpa, inty_dpa, bathyt, &
1329 dz_neglect, usemasswghtinterp)
1332 call mom_error(fatal,
"No analytic integration option is available with this EOS!")
1341 if (.not.
associated(eos))
call mom_error(fatal, &
1342 "query_compressible called with an unassociated EOS_type EOS.")
1348 subroutine eos_init(param_file, EOS, US)
1349 type(param_file_type),
intent(in) :: param_file
1354 #include "version_variable.h" 1355 character(len=40) :: mdl =
"MOM_EOS" 1356 character(len=40) :: tmpstr
1361 call log_version(param_file, mdl, version,
"")
1363 call get_param(param_file, mdl,
"EQN_OF_STATE", tmpstr, &
1364 "EQN_OF_STATE determines which ocean equation of state "//&
1365 "should be used. Currently, the valid choices are "//&
1366 '"LINEAR", "UNESCO", "WRIGHT", "NEMO" and "TEOS10". '//&
1367 "This is only used if USE_EOS is true.", default=
eos_default)
1368 select case (uppercase(tmpstr))
1380 call mom_error(fatal,
"interpret_eos_selection: EQN_OF_STATE "//&
1381 trim(tmpstr) //
"in input file is invalid.")
1383 call mom_mesg(
'interpret_eos_selection: equation of state set to "' // &
1384 trim(tmpstr)//
'"', 5)
1387 eos%Compressible = .false.
1388 call get_param(param_file, mdl,
"RHO_T0_S0", eos%Rho_T0_S0, &
1390 "this is the density at T=0, S=0.", units=
"kg m-3", &
1392 call get_param(param_file, mdl,
"DRHO_DT", eos%dRho_dT, &
1394 "this is the partial derivative of density with "//&
1395 "temperature.", units=
"kg m-3 K-1", default=-0.2)
1396 call get_param(param_file, mdl,
"DRHO_DS", eos%dRho_dS, &
1398 "this is the partial derivative of density with "//&
1399 "salinity.", units=
"kg m-3 PSU-1", default=0.8)
1402 call get_param(param_file, mdl,
"EOS_QUADRATURE", eos%EOS_quadrature, &
1403 "If true, always use the generic (quadrature) code "//&
1404 "code for the integrals of density.", default=.false.)
1406 call get_param(param_file, mdl,
"TFREEZE_FORM", tmpstr, &
1407 "TFREEZE_FORM determines which expression should be "//&
1408 "used for the freezing point. Currently, the valid "//&
1409 'choices are "LINEAR", "MILLERO_78", "TEOS10"', &
1411 select case (uppercase(tmpstr))
1419 call mom_error(fatal,
"interpret_eos_selection: TFREEZE_FORM "//&
1420 trim(tmpstr) //
"in input file is invalid.")
1424 call get_param(param_file, mdl,
"TFREEZE_S0_P0",eos%TFr_S0_P0, &
1426 "this is the freezing potential temperature at "//&
1427 "S=0, P=0.", units=
"deg C", default=0.0)
1428 call get_param(param_file, mdl,
"DTFREEZE_DS",eos%dTFr_dS, &
1430 "this is the derivative of the freezing potential "//&
1431 "temperature with salinity.", &
1432 units=
"deg C PSU-1", default=-0.054)
1433 call get_param(param_file, mdl,
"DTFREEZE_DP",eos%dTFr_dP, &
1435 "this is the derivative of the freezing potential "//&
1436 "temperature with pressure.", &
1437 units=
"deg C Pa-1", default=0.0)
1442 call mom_error(fatal,
"interpret_eos_selection: EOS_TEOS10 or EOS_NEMO \n" //&
1443 "should only be used along with TFREEZE_FORM = TFREEZE_TEOS10 .")
1447 eos%m_to_Z = 1. ;
if (
present(us)) eos%m_to_Z = us%m_to_Z
1448 eos%kg_m3_to_R = 1. ;
if (
present(us)) eos%kg_m3_to_R = us%kg_m3_to_R
1449 eos%R_to_kg_m3 = 1. ;
if (
present(us)) eos%R_to_kg_m3 = us%R_to_kg_m3
1450 eos%RL2_T2_to_Pa = 1. ;
if (
present(us)) eos%RL2_T2_to_Pa = us%RL2_T2_to_Pa
1451 eos%L_T_to_m_s = 1. ;
if (
present(us)) eos%L_T_to_m_s = us%L_T_to_m_s
1456 subroutine eos_manual_init(EOS, form_of_EOS, form_of_TFreeze, EOS_quadrature, Compressible, &
1457 Rho_T0_S0, drho_dT, dRho_dS, TFr_S0_P0, dTFr_dS, dTFr_dp)
1459 integer,
optional,
intent(in) :: form_of_EOS
1460 integer,
optional,
intent(in) :: form_of_TFreeze
1462 logical,
optional,
intent(in) :: EOS_quadrature
1464 logical,
optional,
intent(in) :: Compressible
1465 real ,
optional,
intent(in) :: Rho_T0_S0
1466 real ,
optional,
intent(in) :: drho_dT
1468 real ,
optional,
intent(in) :: dRho_dS
1470 real ,
optional,
intent(in) :: TFr_S0_P0
1471 real ,
optional,
intent(in) :: dTFr_dS
1473 real ,
optional,
intent(in) :: dTFr_dp
1476 if (
present(form_of_eos )) eos%form_of_EOS = form_of_eos
1477 if (
present(form_of_tfreeze)) eos%form_of_TFreeze = form_of_tfreeze
1478 if (
present(eos_quadrature )) eos%EOS_quadrature = eos_quadrature
1479 if (
present(compressible )) eos%Compressible = compressible
1480 if (
present(rho_t0_s0 )) eos%Rho_T0_S0 = rho_t0_s0
1481 if (
present(drho_dt )) eos%drho_dT = drho_dt
1482 if (
present(drho_ds )) eos%dRho_dS = drho_ds
1483 if (
present(tfr_s0_p0 )) eos%TFr_S0_P0 = tfr_s0_p0
1484 if (
present(dtfr_ds )) eos%dTFr_dS = dtfr_ds
1485 if (
present(dtfr_dp )) eos%dTFr_dp = dtfr_dp
1493 if (.not.
associated(eos))
allocate(eos)
1500 if (
associated(eos))
deallocate(eos)
1508 subroutine eos_use_linear(Rho_T0_S0, dRho_dT, dRho_dS, EOS, use_quadrature)
1509 real,
intent(in) :: Rho_T0_S0
1510 real,
intent(in) :: dRho_dT
1511 real,
intent(in) :: dRho_dS
1512 logical,
optional,
intent(in) :: use_quadrature
1516 if (.not.
associated(eos))
call mom_error(fatal, &
1517 "MOM_EOS.F90: EOS_use_linear() called with an unassociated EOS_type EOS.")
1520 eos%Compressible = .false.
1521 eos%Rho_T0_S0 = rho_t0_s0
1522 eos%dRho_dT = drho_dt
1523 eos%dRho_dS = drho_ds
1524 eos%EOS_quadrature = .false.
1525 if (
present(use_quadrature)) eos%EOS_quadrature = use_quadrature
1532 integer,
intent(in) :: kd
1533 type(hor_index_type),
intent(in) :: HI
1534 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed,kd), &
1536 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed,kd), &
1538 real,
dimension(HI%isd:HI%ied,HI%jsd:HI%jed,kd), &
1539 intent(in) :: mask_z
1543 real :: gsw_sr_from_sp, gsw_ct_from_pt, gsw_sa_from_sp
1546 if (.not.
associated(eos))
call mom_error(fatal, &
1547 "convert_temp_salt_to_TEOS10 called with an unassociated EOS_type EOS.")
1551 do k=1,kd ;
do j=hi%jsc,hi%jec ;
do i=hi%isc,hi%iec
1552 if (mask_z(i,j,k) >= 1.0)
then 1553 s(i,j,k) = gsw_sr_from_sp(s(i,j,k))
1558 t(i,j,k) = gsw_ct_from_pt(s(i,j,k), t(i,j,k))
1560 enddo ;
enddo ;
enddo 1572 subroutine extract_member_eos(EOS, form_of_EOS, form_of_TFreeze, EOS_quadrature, Compressible, &
1573 Rho_T0_S0, drho_dT, dRho_dS, TFr_S0_P0, dTFr_dS, dTFr_dp)
1575 integer,
optional,
intent(out) :: form_of_EOS
1576 integer,
optional,
intent(out) :: form_of_TFreeze
1578 logical,
optional,
intent(out) :: EOS_quadrature
1580 logical,
optional,
intent(out) :: Compressible
1581 real ,
optional,
intent(out) :: Rho_T0_S0
1582 real ,
optional,
intent(out) :: drho_dT
1584 real ,
optional,
intent(out) :: dRho_dS
1586 real ,
optional,
intent(out) :: TFr_S0_P0
1587 real ,
optional,
intent(out) :: dTFr_dS
1589 real ,
optional,
intent(out) :: dTFr_dp
1592 if (
present(form_of_eos )) form_of_eos = eos%form_of_EOS
1593 if (
present(form_of_tfreeze)) form_of_tfreeze = eos%form_of_TFreeze
1594 if (
present(eos_quadrature )) eos_quadrature = eos%EOS_quadrature
1595 if (
present(compressible )) compressible = eos%Compressible
1596 if (
present(rho_t0_s0 )) rho_t0_s0 = eos%Rho_T0_S0
1597 if (
present(drho_dt )) drho_dt = eos%drho_dT
1598 if (
present(drho_ds )) drho_ds = eos%dRho_dS
1599 if (
present(tfr_s0_p0 )) tfr_s0_p0 = eos%TFr_S0_P0
1600 if (
present(dtfr_ds )) dtfr_ds = eos%dTFr_dS
1601 if (
present(dtfr_dp )) dtfr_dp = eos%dTFr_dp
subroutine calculate_tfreeze_scalar(S, pressure, T_fr, EOS, pres_scale)
Calls the appropriate subroutine to calculate the freezing point for scalar inputs.
subroutine calculate_density_scalar(T, S, pressure, rho, EOS, rho_ref, scale)
Calls the appropriate subroutine to calculate density of sea water for scalar inputs. If rho_ref is present, the anomaly with respect to rho_ref is returned. The pressure and density can be rescaled with the US. If both the US and scale arguments are present the density scaling uses the product of the two scaling factors.
Calculate the derivatives of specific volume with temperature and salinity from T, S, and P.
integer function, dimension(2), public eos_domain(HI, halo)
This subroutine returns a two point integer array indicating the domain of i-indices to work on in EO...
logical function, public query_compressible(EOS)
Returns true if the equation of state is compressible (i.e. has pressure dependence) ...
subroutine calculate_stanley_density_scalar(T, S, pressure, Tvar, TScov, Svar, rho, EOS, rho_ref, scale)
Calls the appropriate subroutine to calculate density of sea water for scalar inputs including the va...
subroutine calculate_density_array(T, S, pressure, rho, start, npts, EOS, rho_ref, scale)
Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs...
subroutine, public eos_allocate(EOS)
Allocates EOS_type.
Calculates density of sea water from T, S and P.
Calculates the freezing point of sea water from T, S and P.
subroutine, public eos_end(EOS)
Deallocates EOS_type.
subroutine, public convert_temp_salt_for_teos10(T, S, HI, kd, mask_z, EOS)
Convert T&S to Absolute Salinity and Conservative Temperature if using TEOS10.
subroutine calculate_density_derivs_scalar(T, S, pressure, drho_dT, drho_dS, EOS, scale)
Calls the appropriate subroutines to calculate density derivatives by promoting a scalar to a one-ele...
character *(10), parameter eos_unesco_string
A string for specifying the equation of state.
integer, parameter, public eos_unesco
A named integer specifying an equation of state.
integer, parameter, public eos_teos10
A named integer specifying an equation of state.
integer, parameter, public eos_wright
A named integer specifying an equation of state.
subroutine calculate_density_derivs_1d(T, S, pressure, drho_dT, drho_dS, EOS, dom, scale)
Calls the appropriate subroutine to calculate density derivatives for 1-D array inputs.
subroutine calc_spec_vol_derivs_1d(T, S, pressure, dSV_dT, dSV_dS, EOS, dom, scale)
Calls the appropriate subroutine to calculate specific volume derivatives for 1-d array inputs...
character *(10), parameter eos_wright_string
A string for specifying the equation of state.
integer, parameter tfreeze_teos10
A named integer specifying a freezing point expression.
subroutine calculate_spec_vol_array(T, S, pressure, specvol, start, npts, EOS, spv_ref, scale)
Calls the appropriate subroutine to calculate the specific volume of sea water for 1-D array inputs...
character *(10), parameter eos_linear_string
A string for specifying the equation of state.
subroutine calculate_tfreeze_array(S, pressure, T_fr, start, npts, EOS, pres_scale)
Calls the appropriate subroutine to calculate the freezing point for a 1-D array. ...
subroutine calc_spec_vol_1d(T, S, pressure, specvol, EOS, dom, spv_ref, scale)
Calls the appropriate subroutine to calculate the specific volume of sea water for 1-D array inputs...
subroutine calculate_density_second_derivs_scalar(T, S, pressure, drho_dS_dS, drho_dS_dT, drho_dT_dT, drho_dS_dP, drho_dT_dP, EOS, scale)
Calls the appropriate subroutine to calculate density second derivatives for scalar nputs...
subroutine calculate_compress_array(T, S, press, rho, drho_dp, start, npts, EOS)
Calls the appropriate subroutine to calculate the density and compressibility for 1-D array inputs...
Describes various unit conversion factors.
Calculates specific volume of sea water from T, S and P.
Calculates the compressibility of water from T, S, and P.
integer, parameter tfreeze_linear
A named integer specifying a freezing point expression.
integer, parameter tfreeze_millero
A named integer specifying a freezing point expression.
character *(10), parameter tfreeze_teos10_string
A string for specifying the freezing point expression.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
logical function, public eos_quadrature(EOS)
Return value of EOS_quadrature.
subroutine calculate_spec_vol_derivs_array(T, S, pressure, dSV_dT, dSV_dS, start, npts, EOS)
Calls the appropriate subroutine to calculate specific volume derivatives for an array.
subroutine calc_spec_vol_scalar(T, S, pressure, specvol, EOS, spv_ref, scale)
Calls the appropriate subroutine to calculate specific volume of sea water for scalar inputs...
subroutine, public analytic_int_specific_vol_dp(T, S, p_t, p_b, alpha_ref, HI, EOS, dza, intp_dza, intx_dza, inty_dza, halo_size, bathyP, dP_tiny, useMassWghtInterp)
Calls the appropriate subroutine to calculate analytical and nearly-analytical integrals in pressure ...
Calculate the derivatives of density with temperature and salinity from T, S, and P...
Provides subroutines for quantities specific to the equation of state.
integer, parameter, public eos_nemo
A named integer specifying an equation of state.
subroutine calculate_stanley_density_1d(T, S, pressure, Tvar, TScov, Svar, rho, EOS, dom, rho_ref, scale)
Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs including...
subroutine, public eos_init(param_file, EOS, US)
Initializes EOS_type by allocating and reading parameters.
character *(10), parameter tfreeze_millero_string
A string for specifying freezing point expression.
subroutine calculate_compress_scalar(T, S, pressure, rho, drho_dp, EOS)
Calculate density and compressibility for a scalar. This just promotes the scalar to an array with a ...
subroutine calculate_stanley_density_array(T, S, pressure, Tvar, TScov, Svar, rho, start, npts, EOS, rho_ref, scale)
Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs including...
character *(10), parameter tfreeze_default
The default freezing point expression.
subroutine calculate_density_1d(T, S, pressure, rho, EOS, dom, rho_ref, scale)
Calls the appropriate subroutine to calculate the density of sea water for 1-D array inputs...
integer, parameter, public eos_linear
A named integer specifying an equation of state.
subroutine, public eos_use_linear(Rho_T0_S0, dRho_dT, dRho_dS, EOS, use_quadrature)
Set equation of state structure (EOS) to linear with given coefficients.
character *(10), parameter eos_nemo_string
A string for specifying the equation of state.
character *(10), parameter tfreeze_linear_string
A string for specifying the freezing point expression.
subroutine, public eos_manual_init(EOS, form_of_EOS, form_of_TFreeze, EOS_quadrature, Compressible, Rho_T0_S0, drho_dT, dRho_dS, TFr_S0_P0, dTFr_dS, dTFr_dp)
Manually initialized an EOS type (intended for unit testing of routines which need a specific EOS) ...
subroutine calculate_density_second_derivs_array(T, S, pressure, drho_dS_dS, drho_dS_dT, drho_dT_dT, drho_dS_dP, drho_dT_dP, start, npts, EOS, scale)
Calls the appropriate subroutine to calculate density second derivatives for 1-D array inputs...
Calculates the second derivatives of density with various combinations of temperature, salinity, and pressure from T, S and P.
character *(10), parameter eos_teos10_string
A string for specifying the equation of state.
subroutine, public extract_member_eos(EOS, form_of_EOS, form_of_TFreeze, EOS_quadrature, Compressible, Rho_T0_S0, drho_dT, dRho_dS, TFr_S0_P0, dTFr_dS, dTFr_dp)
Extractor routine for the EOS type if the members need to be accessed outside this module...
character *(10), parameter eos_default
The default equation of state.
subroutine, public analytic_int_density_dz(T, S, z_t, z_b, rho_ref, rho_0, G_e, HI, EOS, dpa, intz_dpa, intx_dpa, inty_dpa, bathyT, dz_neglect, useMassWghtInterp)
This subroutine calculates analytical and nearly-analytical integrals of pressure anomalies across la...
subroutine calculate_density_derivs_array(T, S, pressure, drho_dT, drho_dS, start, npts, EOS, scale)
Calls the appropriate subroutine to calculate density derivatives for 1-D array inputs.
A control structure for the equation of state.