FV3 Bundle
ufo_tpsp2ti_mod.F90
Go to the documentation of this file.
2  implicit none
3  private
5 contains
6  subroutine insitu_t_nl(temp_i, temp_p, salt_p, lono, lato, deptho)
7  !==========================================================================
8  !
9  ! Calculates insitu temperature from potential temp and practical salinity
10  !
11  ! Input:
12  ! ------
13  ! temp_p : Potential teperature [deg C]
14  ! salt_p : Practical salinity [psu]
15  ! layer : Layer thickness [m]
16  ! lono : Ref. Conservative Temperature longitude of obs
17  ! lato : Reference sea-surface height latitude "
18  ! deptho : Depth of t,s depth "
19  !
20  ! Output:
21  ! -------
22  ! temp_i : Insitu temperature interpolated at obs depth location [C]
23  ! ti : Insitu temperature at model levels [C]
24  !--------------------------------------------------------------------------
25  !
26  ! z1=0 ----------- Surface
27  ! t1,s1 h1 (thickness)
28  ! z2 -----------
29  ! t2,s2 h2
30  ! z3 -----------
31  ! t3,s3
32  ! .
33  ! .
34  ! .
35  ! zN-1 -----------
36  ! tN-1,sN-1 hN-1
37  ! zN -----------
38  ! tN,sN hN
39  ! zN+1 ----------- Bottom
40  ! ///////////
41 
42  use gsw_mod_kinds
44  use vert_interp_mod
45  use kinds
46 
47  implicit none
48 
49  real(kind=kind_real), intent(in) :: temp_p !< Potential temperature at observation location [C]
50  real(kind=kind_real), intent(in) :: salt_p !< Practical Salinity at observation location [ppt]
51  real(kind=kind_real), intent(in) :: lono, lato, deptho !< Observation location
52  real(kind=kind_real), intent(inout) :: temp_i !< Vertically interpolated Insitu temperature at
53 
54  real(kind_real) :: tp, sp, prso
55 
56 
57  !< Pressure from depth
58  prso = p_from_z(deptho,lato)
59 
60  ! Insitu temperature
61  temp_i = t_from_pt(temp_p,salt_p,prso,lono,lato)
62 
63  if (isnan(temp_i)) temp_i = 0.0
64 
65  end subroutine insitu_t_nl
66 
67  subroutine insitu_t_tl(dtemp_i, dtemp_p, dsalt_p, temp_p, salt_p, lono, lato, deptho, Jacobian)
68 
69  use gsw_mod_kinds
71  use vert_interp_mod
72  use kinds
73 
74  implicit none
75 
76  real(kind=kind_real), intent(in) :: dtemp_p !< Potential temperature at observation location [C]
77  real(kind=kind_real), intent(in) :: dsalt_p !< Practical Salinity at observation location [ppt]
78  real(kind=kind_real), intent(in) :: temp_p !< Bkg potential temperature at observation location [C]
79  real(kind=kind_real), intent(in) :: salt_p !< Bkg practical Salinity at observation location [ppt]
80  real(kind=kind_real), intent(in) :: lono, lato, deptho !< Observation location
81  real(kind=kind_real), intent(inout) :: dtemp_i !< Vertically interpolated Insitu temperature at
82  real(kind=kind_real), intent(in), optional :: jacobian(2) !< Precomputed Jacobian
83 
84  real(kind=kind_real) :: prso, jac(2)
85 
86  !< Pressure from depth
87  prso = p_from_z(deptho,lato)
88 
89  ! Jacobian: dTi/dTp and dTi/dSp
90  if (present(jacobian)) then
91  jac=jacobian
92  else
93  call insitu_t_jac(jac, temp_p, salt_p, lono, lato, deptho)
94  end if
95 
96  ! Tangent of insitu temperature at (tp,sp)
97  dtemp_i = jac(1)*dtemp_p + jac(2)*dsalt_p
98 
99  if (isnan(dtemp_i)) dtemp_i = 0.0
100 
101  end subroutine insitu_t_tl
102 
103  subroutine insitu_t_tlad(dtemp_i, dtemp_p, dsalt_p, temp_p, salt_p, lono, lato, deptho, Jacobian)
104 
105  use gsw_mod_kinds
107  use vert_interp_mod
108  use kinds
109 
110  implicit none
111 
112  real(kind=kind_real), intent(inout) :: dtemp_p !< Potential temperature at observation location [C]
113  real(kind=kind_real), intent(inout) :: dsalt_p !< Practical Salinity at observation location [ppt]
114  real(kind=kind_real), intent(in) :: temp_p !< Bkg potential temperature at observation location [C]
115  real(kind=kind_real), intent(in) :: salt_p !< Bkg practical Salinity at observation location [ppt]
116  real(kind=kind_real), intent(in) :: lono, lato, deptho !< Observation location
117  real(kind=kind_real), intent(in) :: dtemp_i !< Insitu temperature increment at obs loc
118  real(kind=kind_real), intent(in), optional :: jacobian(2) !< Precomputed Jacobian
119 
120  real(kind=kind_real) :: prso, jac(2)
121 
122  !< Pressure from depth
123  prso = p_from_z(deptho,lato)
124 
125  ! Jacobian: dTi/dTp and dTi/dSp
126  if (present(jacobian)) then
127  jac=jacobian
128  else
129  call insitu_t_jac(jac, temp_p, salt_p, lono, lato, deptho)
130  end if
131 
132  !< Adjoint
133  dtemp_p = dtemp_p + jac(1)*dtemp_i
134  dsalt_p = dsalt_p + jac(2)*dtemp_i
135 
136  if (isnan(dtemp_p+dsalt_p)) then
137  dtemp_p = 0.0
138  dsalt_p = 0.0
139  end if
140 
141  end subroutine insitu_t_tlad
142 
143  subroutine insitu_t_jac(jac, temp_p, salt_p, lono, lato, deptho)
144  !==========================================================================
145  ! return jacobian at model levels
146  use gsw_mod_kinds
148  use vert_interp_mod
149  use kinds
150 
151  implicit none
152 
153  real(kind=kind_real), intent(in) :: temp_p !< Potential temperature at observation location [C]
154  real(kind=kind_real), intent(in) :: salt_p !< Practical Salinity at observation location [ppt]
155  real(kind=kind_real), intent(in) :: lono, lato, deptho !< Observation location
156  real(kind=kind_real), intent(out) :: jac(2) !< Jacobian (dti/dtp, dti,dsp)
157 
158  real(kind=kind_real) :: pressure
159  real(kind=kind_real) :: delta=1.0e-10
160  real(kind=kind_real) :: delta_tp, delta_sp
161 
162  ! Vertical interpolation
163  real(kind_real) :: wf
164  integer :: wi
165 
166  delta_tp=delta
167  delta_sp=delta
168 
169  !Pressure from depth
170  pressure=p_from_z(deptho,lato)
171 
172  !Jacobian of Insitu temperature
173  !dti/dtp
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
176  !dti/dsp
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
179 
180  !if (isnan(sum(Jac))) Jac = 0.0
181 
182  end subroutine insitu_t_jac
183 
184 
185 end module ufo_tpsp2ti_mod
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.
Definition: vert_interp.F90:8
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)