6 subroutine insitu_t_nl(temp_i, temp_p, salt_p, lono, lato, deptho)
49 real(kind=kind_real),
intent(in) :: temp_p
50 real(kind=kind_real),
intent(in) :: salt_p
51 real(kind=kind_real),
intent(in) :: lono, lato, deptho
52 real(kind=kind_real),
intent(inout) :: temp_i
54 real(kind_real) :: tp, sp, prso
61 temp_i =
t_from_pt(temp_p,salt_p,prso,lono,lato)
63 if (isnan(temp_i)) temp_i = 0.0
67 subroutine insitu_t_tl(dtemp_i, dtemp_p, dsalt_p, temp_p, salt_p, lono, lato, deptho, Jacobian)
76 real(kind=kind_real),
intent(in) :: dtemp_p
77 real(kind=kind_real),
intent(in) :: dsalt_p
78 real(kind=kind_real),
intent(in) :: temp_p
79 real(kind=kind_real),
intent(in) :: salt_p
80 real(kind=kind_real),
intent(in) :: lono, lato, deptho
81 real(kind=kind_real),
intent(inout) :: dtemp_i
82 real(kind=kind_real),
intent(in),
optional :: jacobian(2)
84 real(kind=kind_real) :: prso, jac(2)
90 if (
present(jacobian))
then 93 call insitu_t_jac(jac, temp_p, salt_p, lono, lato, deptho)
97 dtemp_i = jac(1)*dtemp_p + jac(2)*dsalt_p
99 if (isnan(dtemp_i)) dtemp_i = 0.0
103 subroutine insitu_t_tlad(dtemp_i, dtemp_p, dsalt_p, temp_p, salt_p, lono, lato, deptho, Jacobian)
112 real(kind=kind_real),
intent(inout) :: dtemp_p
113 real(kind=kind_real),
intent(inout) :: dsalt_p
114 real(kind=kind_real),
intent(in) :: temp_p
115 real(kind=kind_real),
intent(in) :: salt_p
116 real(kind=kind_real),
intent(in) :: lono, lato, deptho
117 real(kind=kind_real),
intent(in) :: dtemp_i
118 real(kind=kind_real),
intent(in),
optional :: jacobian(2)
120 real(kind=kind_real) :: prso, jac(2)
126 if (
present(jacobian))
then 129 call insitu_t_jac(jac, temp_p, salt_p, lono, lato, deptho)
133 dtemp_p = dtemp_p + jac(1)*dtemp_i
134 dsalt_p = dsalt_p + jac(2)*dtemp_i
136 if (isnan(dtemp_p+dsalt_p))
then 143 subroutine insitu_t_jac(jac, temp_p, salt_p, lono, lato, deptho)
153 real(kind=kind_real),
intent(in) :: temp_p
154 real(kind=kind_real),
intent(in) :: salt_p
155 real(kind=kind_real),
intent(in) :: lono, lato, deptho
156 real(kind=kind_real),
intent(out) :: jac(2)
158 real(kind=kind_real) :: pressure
159 real(kind=kind_real) :: delta=1.0e-10
160 real(kind=kind_real) :: delta_tp, delta_sp
163 real(kind_real) :: wf
174 jac(1) = (
t_from_pt(temp_p+delta_tp,salt_p,pressure,lono,lato) -&
175 &
t_from_pt(temp_p,salt_p,pressure,lono,lato) )/delta_tp
177 jac(2) = (
t_from_pt(temp_p,salt_p+delta_sp,pressure,lono,lato) -&
178 &
t_from_pt(temp_p,salt_p,pressure,lono,lato) )/delta_tp
subroutine, public insitu_t_jac(jac, temp_p, salt_p, lono, lato, deptho)
subroutine, public insitu_t_tl(dtemp_i, dtemp_p, dsalt_p, temp_p, salt_p, lono, lato, deptho, Jacobian)
subroutine, public insitu_t_tlad(dtemp_i, dtemp_p, dsalt_p, temp_p, salt_p, lono, lato, deptho, Jacobian)
Fortran module to perform linear interpolation.
real(kind_real) function, public t_from_pt(pt_in, sp_in, p_in, lon_in, lat_in)
subroutine, public insitu_t_nl(temp_i, temp_p, salt_p, lono, lato, deptho)
pure real(kind_real) function, public p_from_z(dpth, xlat)