15 implicit none ;
private
17 #include <MOM_memory.h>
19 public user_change_diff, user_change_diff_init, user_change_diff_end
35 logical :: use_abs_lat
47 subroutine user_change_diff(h, tv, G, GV, US, CS, Kd_lay, Kd_int, T_f, S_f, Kd_int_add)
50 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
56 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
optional,
intent(inout) :: kd_lay
58 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
optional,
intent(inout) :: kd_int
60 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
optional,
intent(in) :: t_f
62 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
optional,
intent(in) :: s_f
64 real,
dimension(:,:,:),
optional,
pointer :: kd_int_add
68 real :: rcv(szi_(g),szk_(g))
69 real :: p_ref(szi_(g))
74 logical :: store_kd_add
75 integer,
dimension(2) :: eosdom
76 integer :: i, j, k, is, ie, js, je, nz
77 integer :: isd, ied, jsd, jed
81 character(len=200) :: mesg
83 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
84 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
86 if (.not.
associated(cs))
call mom_error(fatal,
"user_set_diffusivity: "//&
87 "Module must be initialized before it is used.")
89 use_eos =
associated(tv%eqn_of_state)
90 if (.not.use_eos)
return
91 store_kd_add = .false.
92 if (
present(kd_int_add)) store_kd_add =
associated(kd_int_add)
94 if (.not.range_ok(cs%lat_range))
then
95 write(mesg,
'(4(1pe15.6))') cs%lat_range(1:4)
96 call mom_error(fatal,
"user_set_diffusivity: bad latitude range: \n "//&
99 if (.not.range_ok(cs%rho_range))
then
100 write(mesg,
'(4(1pe15.6))') cs%rho_range(1:4)
101 call mom_error(fatal,
"user_set_diffusivity: bad density range: \n "//&
105 if (store_kd_add) kd_int_add(:,:,:) = 0.0
107 do i=is,ie ; p_ref(i) = tv%P_Ref ;
enddo
108 eosdom(:) = eos_domain(g%HI)
110 if (
present(t_f) .and.
present(s_f))
then
112 call calculate_density(t_f(:,j,k), s_f(:,j,k), p_ref, rcv(:,k), tv%eqn_of_state, eosdom)
116 call calculate_density(tv%T(:,j,k), tv%S(:,j,k), p_ref, rcv(:,k), tv%eqn_of_state, eosdom)
120 if (
present(kd_lay))
then
121 do k=1,nz ;
do i=is,ie
122 if (cs%use_abs_lat)
then
123 lat_fn = val_weights(abs(g%geoLatT(i,j)), cs%lat_range)
125 lat_fn = val_weights(g%geoLatT(i,j), cs%lat_range)
127 rho_fn = val_weights(rcv(i,k), cs%rho_range)
128 if (rho_fn * lat_fn > 0.0) &
129 kd_lay(i,j,k) = kd_lay(i,j,k) + cs%Kd_add * rho_fn * lat_fn
132 if (
present(kd_int))
then
133 do k=2,nz ;
do i=is,ie
134 if (cs%use_abs_lat)
then
135 lat_fn = val_weights(abs(g%geoLatT(i,j)), cs%lat_range)
137 lat_fn = val_weights(g%geoLatT(i,j), cs%lat_range)
139 rho_fn = val_weights( 0.5*(rcv(i,k-1) + rcv(i,k)), cs%rho_range)
140 if (rho_fn * lat_fn > 0.0)
then
141 kd_int(i,j,k) = kd_int(i,j,k) + cs%Kd_add * rho_fn * lat_fn
142 if (store_kd_add) kd_int_add(i,j,k) = cs%Kd_add * rho_fn * lat_fn
148 end subroutine user_change_diff
151 function range_ok(range)
result(OK)
152 real,
dimension(4),
intent(in) :: range
155 ok = ((range(1) <= range(2)) .and. (range(2) <= range(3)) .and. &
156 (range(3) <= range(4)))
158 end function range_ok
165 function val_weights(val, range)
result(ans)
166 real,
intent(in) :: val
167 real,
dimension(4),
intent(in) :: range
173 if ((val > range(1)) .and. (val < range(4)))
then
174 if (val < range(2))
then
176 x = (val - range(1)) / (range(2) - range(1))
177 ans = x**2 * (3.0 - 2.0 * x)
178 elseif (val > range(3))
then
180 x = (range(4) - val) / (range(4) - range(3))
181 ans = x**2 * (3.0 - 2.0 * x)
187 end function val_weights
190 subroutine user_change_diff_init(Time, G, GV, US, param_file, diag, CS)
191 type(time_type),
intent(in) :: time
198 type(
diag_ctrl),
target,
intent(inout) :: diag
205 #include "version_variable.h"
206 character(len=40) :: mdl =
"user_set_diffusivity"
207 character(len=200) :: mesg
208 integer :: i, j, is, ie, js, je
210 if (
associated(cs))
then
211 call mom_error(warning,
"diabatic_entrain_init called with an associated "// &
212 "control structure.")
217 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
223 call get_param(param_file, mdl,
"USER_KD_ADD", cs%Kd_add, &
224 "A user-specified additional diffusivity over a range of "//&
225 "latitude and density.", default=0.0, units=
"m2 s-1", &
226 scale=us%m2_s_to_Z2_T)
227 if (cs%Kd_add /= 0.0)
then
228 call get_param(param_file, mdl,
"USER_KD_ADD_LAT_RANGE", cs%lat_range(:), &
229 "Four successive values that define a range of latitudes "//&
230 "over which the user-specified extra diffusivity is "//&
231 "applied. The four values specify the latitudes at "//&
232 "which the extra diffusivity starts to increase from 0, "//&
233 "hits its full value, starts to decrease again, and is "//&
234 "back to 0.", units=
"degree", default=-1.0e9)
235 call get_param(param_file, mdl,
"USER_KD_ADD_RHO_RANGE", cs%rho_range(:), &
236 "Four successive values that define a range of potential "//&
237 "densities over which the user-given extra diffusivity "//&
238 "is applied. The four values specify the density at "//&
239 "which the extra diffusivity starts to increase from 0, "//&
240 "hits its full value, starts to decrease again, and is "//&
241 "back to 0.", units=
"kg m-3", default=-1.0e9, scale=us%kg_m3_to_R)
242 call get_param(param_file, mdl,
"USER_KD_ADD_USE_ABS_LAT", cs%use_abs_lat, &
243 "If true, use the absolute value of latitude when "//&
244 "checking whether a point fits into range of latitudes.", &
248 if (.not.range_ok(cs%lat_range))
then
249 write(mesg,
'(4(1pe15.6))') cs%lat_range(1:4)
250 call mom_error(fatal,
"user_set_diffusivity: bad latitude range: \n "//&
253 if (.not.range_ok(cs%rho_range))
then
254 write(mesg,
'(4(1pe15.6))') cs%rho_range(1:4)
255 call mom_error(fatal,
"user_set_diffusivity: bad density range: \n "//&
259 end subroutine user_change_diff_init
262 subroutine user_change_diff_end(CS)
266 if (
associated(cs))
deallocate(cs)
268 end subroutine user_change_diff_end