FV3 Bundle
xbt_drop_rate_adjust.f90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 module xbt_adjust
20 
22 
23  implicit none
24 
25  real, parameter :: s1=0.30731408,s2=6.707e-9,s3=-8.1899e-5,sa=3.227,&
26  sb=-2.17e-4
27  real, parameter :: t1=0.29585798,t2=1.002e-9,t3=-3.1658e-5,ta=3.426,&
28  tb=-4.7e-4
29 
30 contains
31 
32 subroutine xbt_drop_rate_adjust(station)
33 
34 
35  type(ocean_profile_type), intent(inout) :: station
36 
37 
38  integer :: k
39  real :: dpth_orig, dpth_new,tdrop
40  integer :: fix_depth
41 
42  if (.not. station%accepted) return
43 
44  fix_depth = int(station%fix_depth)
45  select case(fix_depth)
46  case(-1)
47  return
48  case(0)
49  return
50  case(1)
51 ! use Hanawa et al (1994) drop rate correction
52  do k=1,station%levels
53  dpth_orig = station%depth(k)
54  if (dpth_orig .ne. missing_value) then
55  dpth_new = (1.0417*dpth_orig) - (75.096*(1.0-((1.0-(0.0002063*dpth_orig)))**0.5))
56  station%depth(k)=dpth_new
57  endif
58  enddo
59  station%fix_depth=-1.0
60  case (2)
61 ! use Kizu et al (2005) correction
62  do k=1,station%levels
63  dpth_orig = station%depth(k)
64  if (dpth_orig .ne. missing_value) then
65  if (dpth_orig .le. 250.0) then
66  dpth_new = dpth_orig*0.9572
67  else if (dpth_orig .le. 500.) then
68  dpth_new = dpth_orig*0.9565
69  else if (dpth_orig .le. 750.0) then
70  dpth_new = dpth_orig*0.9558
71  else if (dpth_orig .le. 1000.) then
72  dpth_new = dpth_orig*0.9550
73  else if (dpth_orig .le. 1250.0) then
74  dpth_new = dpth_orig*0.9542
75  else if (dpth_orig .le. 1500.0) then
76  dpth_new = dpth_orig*0.9533
77  else
78  dpth_new = dpth_orig*0.9524
79  endif
80  station%depth(k)=dpth_new
81  endif
82  enddo
83  station%fix_depth=-1.0
84  case(103)
85  do k=1,station%levels
86  dpth_orig = station%depth(k)
87  if (dpth_orig .ne. missing_value) then
88  tdrop=(s1*dpth_orig + s2) - s3
89  dpth_new = sa*tdrop + sb*tdrop*tdrop
90  station%depth(k)=dpth_new
91  endif
92  enddo
93  station%fix_depth=-1.0
94  case(104)
95  do k=1,station%levels
96  dpth_orig = station%depth(k)
97  if (dpth_orig .ne. missing_value) then
98  tdrop=(t1*dpth_orig + t2) - t3
99  dpth_new = ta*tdrop + tb*tdrop*tdrop
100  station%depth(k)=dpth_new
101  endif
102  enddo
103  station%fix_depth=-1.0
104  end select
105 
106  return
107 end subroutine xbt_drop_rate_adjust
108 
109 end module xbt_adjust
real, parameter s2
real, parameter t3
real, parameter sa
real, parameter ta
real, parameter s1
real, parameter t2
real, parameter sb
real, parameter s3
real, parameter t1
subroutine xbt_drop_rate_adjust(station)
real, parameter tb
real, parameter, public missing_value
Definition: oda_types.F90:69