6 use mom_cpu_clock,
only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
13 use mom_eos,
only : extract_member_eos, eos_linear, eos_teos10, eos_wright
19 use mom_remapping,
only : extract_member_remapping_cs, build_reconstructions_1d
20 use mom_remapping,
only : average_value_ppoly, remappingschemesdoc, remappingdefaultscheme
25 use ppm_functions,
only : ppm_reconstruction, ppm_boundary_extrapolation
32 use iso_fortran_env,
only : stdout=>output_unit, stderr=>error_unit
34 implicit none ;
private
36 #include <MOM_memory.h>
38 public neutral_diffusion, neutral_diffusion_init, neutral_diffusion_end
39 public neutral_diffusion_calc_coeffs
40 public neutral_diffusion_unit_tests
47 logical :: continuous_reconstruction = .true.
48 logical :: debug = .false.
49 logical :: hard_fail_heff
55 logical :: interior_only
58 real,
allocatable,
dimension(:,:,:) :: upol
59 real,
allocatable,
dimension(:,:,:) :: upor
60 integer,
allocatable,
dimension(:,:,:) :: ukol
62 integer,
allocatable,
dimension(:,:,:) :: ukor
64 real,
allocatable,
dimension(:,:,:) :: uheff
65 real,
allocatable,
dimension(:,:,:) :: vpol
66 real,
allocatable,
dimension(:,:,:) :: vpor
67 integer,
allocatable,
dimension(:,:,:) :: vkol
69 integer,
allocatable,
dimension(:,:,:) :: vkor
71 real,
allocatable,
dimension(:,:,:) :: vheff
73 real,
allocatable,
dimension(:,:,:,:) :: ppoly_coeffs_t
74 real,
allocatable,
dimension(:,:,:,:) :: ppoly_coeffs_s
76 real,
allocatable,
dimension(:,:,:) :: drdt
77 real,
allocatable,
dimension(:,:,:) :: drds
78 real,
allocatable,
dimension(:,:,:) :: tint
79 real,
allocatable,
dimension(:,:,:) :: sint
80 real,
allocatable,
dimension(:,:,:) :: pint
82 real,
allocatable,
dimension(:,:,:,:) :: t_i
83 real,
allocatable,
dimension(:,:,:,:) :: s_i
84 real,
allocatable,
dimension(:,:,:,:) :: p_i
85 real,
allocatable,
dimension(:,:,:,:) :: drdt_i
86 real,
allocatable,
dimension(:,:,:,:) :: drds_i
87 integer,
allocatable,
dimension(:,:) :: ns
88 logical,
allocatable,
dimension(:,:,:) :: stable_cell
89 real :: r_to_kg_m3 = 1.0
93 integer :: neutral_pos_method
94 character(len=40) :: delta_rho_form
97 integer :: id_uheff_2d = -1
98 integer :: id_vheff_2d = -1
102 logical :: remap_answers_2018
105 type(
kpp_cs),
pointer :: kpp_csp => null()
110 #include "version_variable.h"
111 character(len=40) :: mdl =
"MOM_neutral_diffusion"
116 logical function neutral_diffusion_init(Time, G, US, param_file, diag, EOS, diabatic_CSp, CS)
117 type(time_type),
target,
intent(in) :: time
120 type(
diag_ctrl),
target,
intent(inout) :: diag
122 type(
eos_type),
target,
intent(in) :: eos
127 character(len=256) :: mesg
128 character(len=80) :: string
129 logical :: default_2018_answers
130 logical :: boundary_extrap
132 if (
associated(cs))
then
133 call mom_error(fatal,
"neutral_diffusion_init called with associated control structure.")
138 call get_param(param_file, mdl,
"USE_NEUTRAL_DIFFUSION", neutral_diffusion_init, &
139 default=.false., do_not_log=.true.)
141 "This module implements neutral diffusion of tracers", &
142 all_default=.not.neutral_diffusion_init)
143 call get_param(param_file, mdl,
"USE_NEUTRAL_DIFFUSION", neutral_diffusion_init, &
144 "If true, enables the neutral diffusion module.", &
147 if (.not.neutral_diffusion_init)
return
155 call get_param(param_file, mdl,
"NDIFF_CONTINUOUS", cs%continuous_reconstruction, &
156 "If true, uses a continuous reconstruction of T and S when "//&
157 "finding neutral surfaces along which diffusion will happen. "//&
158 "If false, a PPM discontinuous reconstruction of T and S "//&
159 "is done which results in a higher order routine but exacts "//&
160 "a higher computational cost.", default=.true.)
161 call get_param(param_file, mdl,
"NDIFF_REF_PRES", cs%ref_pres, &
162 "The reference pressure (Pa) used for the derivatives of "//&
163 "the equation of state. If negative (default), local pressure is used.", &
164 units=
"Pa", default = -1., scale=us%kg_m3_to_R*us%m_s_to_L_T**2)
165 call get_param(param_file, mdl,
"NDIFF_INTERIOR_ONLY", cs%interior_only, &
166 "If true, only applies neutral diffusion in the ocean interior."//&
167 "That is, the algorithm will exclude the surface and bottom"//&
168 "boundary layers.", default = .false.)
171 if ( .not.cs%continuous_reconstruction )
then
172 call get_param(param_file, mdl,
"NDIFF_BOUNDARY_EXTRAP", boundary_extrap, &
173 "Extrapolate at the top and bottommost cells, otherwise \n"// &
174 "assume boundaries are piecewise constant", &
176 call get_param(param_file, mdl,
"NDIFF_REMAPPING_SCHEME", string, &
177 "This sets the reconstruction scheme used "//&
178 "for vertical remapping for all variables. "//&
179 "It can be one of the following schemes: "//&
180 trim(remappingschemesdoc), default=remappingdefaultscheme)
181 call get_param(param_file, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
182 "This sets the default value for the various _2018_ANSWERS parameters.", &
184 call get_param(param_file, mdl,
"REMAPPING_2018_ANSWERS", cs%remap_answers_2018, &
185 "If true, use the order of arithmetic and expressions that recover the "//&
186 "answers from the end of 2018. Otherwise, use updated and more robust "//&
187 "forms of the same expressions.", default=default_2018_answers)
188 call initialize_remapping( cs%remap_CS, string, boundary_extrapolation=boundary_extrap, &
189 answers_2018=cs%remap_answers_2018 )
190 call extract_member_remapping_cs(cs%remap_CS, degree=cs%deg)
191 call get_param(param_file, mdl,
"NEUTRAL_POS_METHOD", cs%neutral_pos_method, &
192 "Method used to find the neutral position \n"// &
193 "1. Delta_rho varies linearly, find 0 crossing \n"// &
194 "2. Alpha and beta vary linearly from top to bottom, \n"// &
195 " Newton's method for neutral position \n"// &
196 "3. Full nonlinear equation of state, use regula falsi \n"// &
197 " for neutral position", default=3)
198 if (cs%neutral_pos_method > 4 .or. cs%neutral_pos_method < 0)
then
199 call mom_error(fatal,
"Invalid option for NEUTRAL_POS_METHOD")
202 call get_param(param_file, mdl,
"DELTA_RHO_FORM", cs%delta_rho_form, &
203 "Determine how the difference in density is calculated \n"// &
204 " full : Difference of in-situ densities \n"// &
205 " no_pressure: Calculated from dRdT, dRdS, but no \n"// &
206 " pressure dependence", &
207 default=
"mid_pressure")
208 if (cs%neutral_pos_method > 1)
then
209 call get_param(param_file, mdl,
"NDIFF_DRHO_TOL", cs%drho_tol, &
210 "Sets the convergence criterion for finding the neutral\n"// &
211 "position within a layer in kg m-3.", &
212 default=1.e-10, scale=us%kg_m3_to_R)
213 call get_param(param_file, mdl,
"NDIFF_X_TOL", cs%x_tol, &
214 "Sets the convergence criterion for a change in nondim\n"// &
215 "position within a layer.", &
217 call get_param(param_file, mdl,
"NDIFF_MAX_ITER", cs%max_iter, &
218 "The maximum number of iterations to be done before \n"// &
219 "exiting the iterative loop to find the neutral surface", &
222 call get_param(param_file, mdl,
"NDIFF_DEBUG", cs%debug, &
223 "Turns on verbose output for discontinuous neutral "//&
224 "diffusion routines.", &
226 call get_param(param_file, mdl,
"HARD_FAIL_HEFF", cs%hard_fail_heff, &
227 "Bring down the model if a problem with heff is detected",&
232 cs%R_to_kg_m3 = us%R_to_kg_m3
234 if (cs%interior_only)
then
235 call extract_diabatic_member(diabatic_csp, kpp_csp=cs%KPP_CSp)
236 call extract_diabatic_member(diabatic_csp, energetic_pbl_csp=cs%energetic_PBL_CSp)
237 if ( .not.
ASSOCIATED(cs%energetic_PBL_CSp) .and. .not.
ASSOCIATED(cs%KPP_CSp) )
then
238 call mom_error(fatal,
"NDIFF_INTERIOR_ONLY is true, but no valid boundary layer scheme was found")
246 if (cs%continuous_reconstruction)
then
248 allocate(cs%dRdT(szi_(g),szj_(g),szk_(g)+1)) ; cs%dRdT(:,:,:) = 0.
249 allocate(cs%dRdS(szi_(g),szj_(g),szk_(g)+1)) ; cs%dRdS(:,:,:) = 0.
252 allocate(cs%T_i(szi_(g),szj_(g),szk_(g),2)) ; cs%T_i(:,:,:,:) = 0.
253 allocate(cs%S_i(szi_(g),szj_(g),szk_(g),2)) ; cs%S_i(:,:,:,:) = 0.
254 allocate(cs%P_i(szi_(g),szj_(g),szk_(g),2)) ; cs%P_i(:,:,:,:) = 0.
255 allocate(cs%dRdT_i(szi_(g),szj_(g),szk_(g),2)) ; cs%dRdT_i(:,:,:,:) = 0.
256 allocate(cs%dRdS_i(szi_(g),szj_(g),szk_(g),2)) ; cs%dRdS_i(:,:,:,:) = 0.
257 allocate(cs%ppoly_coeffs_T(szi_(g),szj_(g),szk_(g),cs%deg+1)) ; cs%ppoly_coeffs_T(:,:,:,:) = 0.
258 allocate(cs%ppoly_coeffs_S(szi_(g),szj_(g),szk_(g),cs%deg+1)) ; cs%ppoly_coeffs_S(:,:,:,:) = 0.
259 allocate(cs%ns(szi_(g),szj_(g))) ; cs%ns(:,:) = 0.
262 allocate(cs%Tint(szi_(g),szj_(g),szk_(g)+1)) ; cs%Tint(:,:,:) = 0.
263 allocate(cs%Sint(szi_(g),szj_(g),szk_(g)+1)) ; cs%Sint(:,:,:) = 0.
264 allocate(cs%Pint(szi_(g),szj_(g),szk_(g)+1)) ; cs%Pint(:,:,:) = 0.
265 allocate(cs%stable_cell(szi_(g),szj_(g),szk_(g))) ; cs%stable_cell(:,:,:) = .true.
267 allocate(cs%uPoL(g%isd:g%ied,g%jsd:g%jed, cs%nsurf)); cs%uPoL(g%isc-1:g%iec,g%jsc:g%jec,:) = 0.
268 allocate(cs%uPoR(g%isd:g%ied,g%jsd:g%jed, cs%nsurf)); cs%uPoR(g%isc-1:g%iec,g%jsc:g%jec,:) = 0.
269 allocate(cs%uKoL(g%isd:g%ied,g%jsd:g%jed, cs%nsurf)); cs%uKoL(g%isc-1:g%iec,g%jsc:g%jec,:) = 0
270 allocate(cs%uKoR(g%isd:g%ied,g%jsd:g%jed, cs%nsurf)); cs%uKoR(g%isc-1:g%iec,g%jsc:g%jec,:) = 0
271 allocate(cs%uHeff(g%isd:g%ied,g%jsd:g%jed,cs%nsurf-1)); cs%uHeff(g%isc-1:g%iec,g%jsc:g%jec,:) = 0
273 allocate(cs%vPoL(g%isd:g%ied,g%jsd:g%jed, cs%nsurf)); cs%vPoL(g%isc:g%iec,g%jsc-1:g%jec,:) = 0.
274 allocate(cs%vPoR(g%isd:g%ied,g%jsd:g%jed, cs%nsurf)); cs%vPoR(g%isc:g%iec,g%jsc-1:g%jec,:) = 0.
275 allocate(cs%vKoL(g%isd:g%ied,g%jsd:g%jed, cs%nsurf)); cs%vKoL(g%isc:g%iec,g%jsc-1:g%jec,:) = 0
276 allocate(cs%vKoR(g%isd:g%ied,g%jsd:g%jed, cs%nsurf)); cs%vKoR(g%isc:g%iec,g%jsc-1:g%jec,:) = 0
277 allocate(cs%vHeff(g%isd:g%ied,g%jsd:g%jed,cs%nsurf-1)); cs%vHeff(g%isc:g%iec,g%jsc-1:g%jec,:) = 0
279 end function neutral_diffusion_init
283 subroutine neutral_diffusion_calc_coeffs(G, GV, US, h, T, S, CS, p_surf)
287 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
288 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: t
289 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: s
291 real,
dimension(SZI_(G),SZJ_(G)),
optional,
intent(in) :: p_surf
295 integer,
dimension(2) :: eosdom
298 real,
dimension(SZK_(G),2) :: ppoly_r_s
299 real,
dimension(SZI_(G), SZJ_(G)) :: heff_sum
300 real,
dimension(SZI_(G),SZJ_(G)) :: hbl
302 real,
dimension(SZI_(G)) :: ref_pres
303 real,
dimension(SZI_(G)) :: rho_tmp
304 real :: h_neglect, h_neglect_edge
305 integer,
dimension(SZI_(G), SZJ_(G)) :: k_top
306 real,
dimension(SZI_(G), SZJ_(G)) :: zeta_top
308 integer,
dimension(SZI_(G), SZJ_(G)) :: k_bot
309 real,
dimension(SZI_(G), SZJ_(G)) :: zeta_bot
312 pa_to_h = 1. / (gv%H_to_RZ * gv%g_Earth)
314 k_top(:,:) = 1 ; k_bot(:,:) = 1
315 zeta_top(:,:) = 0. ; zeta_bot(:,:) = 0.
318 if (cs%interior_only)
then
320 if (
ASSOCIATED(cs%KPP_CSp))
call kpp_get_bld(cs%KPP_CSp, hbl, g, us, m_to_bld_units=gv%m_to_H)
321 if (
ASSOCIATED(cs%energetic_PBL_CSp)) &
322 call energetic_pbl_get_mld(cs%energetic_PBL_CSp, hbl, g, us, m_to_mld_units=gv%m_to_H)
325 do j=g%jsc-1, g%jec+1 ;
do i=g%isc-1,g%iec+1
326 call boundary_k_range(surface, g%ke, h(i,j,:), hbl(i,j), k_top(i,j), zeta_top(i,j), k_bot(i,j), zeta_bot(i,j))
331 if (.not.cs%remap_answers_2018)
then
332 h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
333 elseif (gv%Boussinesq)
then
334 h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
336 h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
340 if (cs%ref_pres>=0.)
then
341 ref_pres(:) = cs%ref_pres
344 if (cs%continuous_reconstruction)
then
350 cs%dRdT_i(:,:,:,:) = 0.
351 cs%dRdS_i(:,:,:,:) = 0.
353 cs%stable_cell(:,:,:) = .true.
357 if (
present(p_surf))
then
358 do j=g%jsc-1,g%jec+1 ;
do i=g%isc-1,g%iec+1
359 cs%Pint(i,j,1) = p_surf(i,j)
364 do k=1,g%ke ;
do j=g%jsc-1,g%jec+1 ;
do i=g%isc-1,g%iec+1
365 cs%Pint(i,j,k+1) = cs%Pint(i,j,k) + h(i,j,k)*(gv%g_Earth*gv%H_to_RZ)
366 enddo ;
enddo ;
enddo
370 if (.not. cs%continuous_reconstruction)
then
371 if (
present(p_surf))
then
372 do j=g%jsc-1,g%jec+1 ;
do i=g%isc-1,g%iec+1
373 cs%P_i(i,j,1,1) = p_surf(i,j)
374 cs%P_i(i,j,1,2) = p_surf(i,j) + h(i,j,1)*(gv%H_to_RZ*gv%g_Earth)
377 do j=g%jsc-1,g%jec+1 ;
do i=g%isc-1,g%iec+1
379 cs%P_i(i,j,1,2) = h(i,j,1)*(gv%H_to_RZ*gv%g_Earth)
382 do k=2,g%ke ;
do j=g%jsc-1,g%jec+1 ;
do i=g%isc-1,g%iec+1
383 cs%P_i(i,j,k,1) = cs%P_i(i,j,k-1,2)
384 cs%P_i(i,j,k,2) = cs%P_i(i,j,k-1,2) + h(i,j,k)*(gv%H_to_RZ*gv%g_Earth)
385 enddo ;
enddo ;
enddo
388 eosdom(:) = eos_domain(g%HI, halo=1)
389 do j = g%jsc-1, g%jec+1
391 do i = g%isc-1, g%iec+1
392 if (cs%continuous_reconstruction)
then
393 call interface_scalar(g%ke, h(i,j,:), t(i,j,:), cs%Tint(i,j,:), 2, h_neglect)
394 call interface_scalar(g%ke, h(i,j,:), s(i,j,:), cs%Sint(i,j,:), 2, h_neglect)
396 call build_reconstructions_1d( cs%remap_CS, g%ke, h(i,j,:), t(i,j,:), cs%ppoly_coeffs_T(i,j,:,:), &
397 cs%T_i(i,j,:,:), ppoly_r_s, imethod, h_neglect, h_neglect_edge )
398 call build_reconstructions_1d( cs%remap_CS, g%ke, h(i,j,:), s(i,j,:), cs%ppoly_coeffs_S(i,j,:,:), &
399 cs%S_i(i,j,:,:), ppoly_r_s, imethod, h_neglect, h_neglect_edge )
403 cs%T_i(i,j,k,1) = evaluation_polynomial( cs%ppoly_coeffs_T(i,j,k,:), cs%deg+1, 0. )
404 cs%T_i(i,j,k,2) = evaluation_polynomial( cs%ppoly_coeffs_T(i,j,k,:), cs%deg+1, 1. )
405 cs%S_i(i,j,k,1) = evaluation_polynomial( cs%ppoly_coeffs_S(i,j,k,:), cs%deg+1, 0. )
406 cs%S_i(i,j,k,2) = evaluation_polynomial( cs%ppoly_coeffs_S(i,j,k,:), cs%deg+1, 1. )
412 if (cs%continuous_reconstruction)
then
414 if (cs%ref_pres<0) ref_pres(:) = cs%Pint(:,j,k)
415 call calculate_density_derivs(cs%Tint(:,j,k), cs%Sint(:,j,k), ref_pres, cs%dRdT(:,j,k), &
416 cs%dRdS(:,j,k), cs%EOS, eosdom)
420 if (cs%ref_pres<0) ref_pres(:) = cs%Pint(:,j,k)
422 call calculate_density_derivs(cs%T_i(:,j,k,1), cs%S_i(:,j,k,1), ref_pres, cs%dRdT_i(:,j,k,1), &
423 cs%dRdS_i(:,j,k,1), cs%EOS, eosdom)
424 if (cs%ref_pres<0) ref_pres(:) = cs%Pint(:,j,k+1)
426 call calculate_density_derivs(cs%T_i(:,j,k,2), cs%S_i(:,j,k,2), ref_pres, cs%dRdT_i(:,j,k,2), &
427 cs%dRdS_i(:,j,k,2), cs%EOS, eosdom)
432 if (.not. cs%continuous_reconstruction)
then
433 do j = g%jsc-1, g%jec+1 ;
do i = g%isc-1, g%iec+1
434 call mark_unstable_cells( cs, g%ke, cs%T_i(i,j,:,:), cs%S_i(i,j,:,:), cs%P_i(i,j,:,:), cs%stable_cell(i,j,:) )
435 if (cs%interior_only)
then
436 if (.not. cs%stable_cell(i,j,k_bot(i,j))) zeta_bot(i,j) = -1.
438 do k = 1, k_bot(i,j)-1
439 cs%stable_cell(i,j,k) = .false.
457 do j = g%jsc, g%jec ;
do i = g%isc-1, g%iec
458 if (g%mask2dCu(i,j) > 0.)
then
459 if (cs%continuous_reconstruction)
then
460 call find_neutral_surface_positions_continuous(g%ke, &
461 cs%Pint(i,j,:), cs%Tint(i,j,:), cs%Sint(i,j,:), cs%dRdT(i,j,:), cs%dRdS(i,j,:), &
462 cs%Pint(i+1,j,:), cs%Tint(i+1,j,:), cs%Sint(i+1,j,:), cs%dRdT(i+1,j,:), cs%dRdS(i+1,j,:), &
463 cs%uPoL(i,j,:), cs%uPoR(i,j,:), cs%uKoL(i,j,:), cs%uKoR(i,j,:), cs%uhEff(i,j,:), &
464 k_bot(i,j), k_bot(i+1,j), zeta_bot(i,j), zeta_bot(i+1,j))
466 call find_neutral_surface_positions_discontinuous(cs, g%ke, &
467 cs%P_i(i,j,:,:), h(i,j,:), cs%T_i(i,j,:,:), cs%S_i(i,j,:,:), cs%ppoly_coeffs_T(i,j,:,:), &
468 cs%ppoly_coeffs_S(i,j,:,:),cs%stable_cell(i,j,:), &
469 cs%P_i(i+1,j,:,:), h(i+1,j,:), cs%T_i(i+1,j,:,:), cs%S_i(i+1,j,:,:), cs%ppoly_coeffs_T(i+1,j,:,:), &
470 cs%ppoly_coeffs_S(i+1,j,:,:), cs%stable_cell(i+1,j,:), &
471 cs%uPoL(i,j,:), cs%uPoR(i,j,:), cs%uKoL(i,j,:), cs%uKoR(i,j,:), cs%uhEff(i,j,:), &
472 hard_fail_heff = cs%hard_fail_heff)
478 do j = g%jsc-1, g%jec ;
do i = g%isc, g%iec
479 if (g%mask2dCv(i,j) > 0.)
then
480 if (cs%continuous_reconstruction)
then
481 call find_neutral_surface_positions_continuous(g%ke, &
482 cs%Pint(i,j,:), cs%Tint(i,j,:), cs%Sint(i,j,:), cs%dRdT(i,j,:), cs%dRdS(i,j,:), &
483 cs%Pint(i,j+1,:), cs%Tint(i,j+1,:), cs%Sint(i,j+1,:), cs%dRdT(i,j+1,:), cs%dRdS(i,j+1,:), &
484 cs%vPoL(i,j,:), cs%vPoR(i,j,:), cs%vKoL(i,j,:), cs%vKoR(i,j,:), cs%vhEff(i,j,:), &
485 k_bot(i,j), k_bot(i,j+1), zeta_bot(i,j), zeta_bot(i,j+1))
487 call find_neutral_surface_positions_discontinuous(cs, g%ke, &
488 cs%P_i(i,j,:,:), h(i,j,:), cs%T_i(i,j,:,:), cs%S_i(i,j,:,:), cs%ppoly_coeffs_T(i,j,:,:), &
489 cs%ppoly_coeffs_S(i,j,:,:),cs%stable_cell(i,j,:), &
490 cs%P_i(i,j+1,:,:), h(i,j+1,:), cs%T_i(i,j+1,:,:), cs%S_i(i,j+1,:,:), cs%ppoly_coeffs_T(i,j+1,:,:), &
491 cs%ppoly_coeffs_S(i,j+1,:,:), cs%stable_cell(i,j+1,:), &
492 cs%vPoL(i,j,:), cs%vPoR(i,j,:), cs%vKoL(i,j,:), cs%vKoR(i,j,:), cs%vhEff(i,j,:), &
493 hard_fail_heff = cs%hard_fail_heff)
502 if (cs%continuous_reconstruction)
then
503 do k = 1, cs%nsurf-1 ;
do j = g%jsc, g%jec ;
do i = g%isc-1, g%iec
504 if (g%mask2dCu(i,j) > 0.) cs%uhEff(i,j,k) = cs%uhEff(i,j,k) * pa_to_h
505 enddo ;
enddo ;
enddo
506 do k = 1, cs%nsurf-1 ;
do j = g%jsc-1, g%jec ;
do i = g%isc, g%iec
507 if (g%mask2dCv(i,j) > 0.) cs%vhEff(i,j,k) = cs%vhEff(i,j,k) * pa_to_h
508 enddo ;
enddo ;
enddo
511 if (cs%id_uhEff_2d>0)
then
513 do k = 1,cs%nsurf-1 ;
do j=g%jsc,g%jec ;
do i=g%isc-1,g%iec
514 heff_sum(i,j) = heff_sum(i,j) + cs%uhEff(i,j,k)
515 enddo ;
enddo ;
enddo
516 call post_data(cs%id_uhEff_2d, heff_sum, cs%diag)
518 if (cs%id_vhEff_2d>0)
then
520 do k = 1,cs%nsurf-1 ;
do j=g%jsc-1,g%jec ;
do i=g%isc,g%iec
521 heff_sum(i,j) = heff_sum(i,j) + cs%vhEff(i,j,k)
522 enddo ;
enddo ;
enddo
523 call post_data(cs%id_vhEff_2d, heff_sum, cs%diag)
526 end subroutine neutral_diffusion_calc_coeffs
529 subroutine neutral_diffusion(G, GV, h, Coef_x, Coef_y, dt, Reg, US, CS)
532 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
533 real,
dimension(SZIB_(G),SZJ_(G)),
intent(in) :: coef_x
534 real,
dimension(SZI_(G),SZJB_(G)),
intent(in) :: coef_y
535 real,
intent(in) :: dt
542 real,
dimension(SZIB_(G),SZJ_(G),CS%nsurf-1) :: uflx
543 real,
dimension(SZI_(G),SZJB_(G),CS%nsurf-1) :: vflx
545 real,
dimension(SZI_(G),SZJ_(G),G%ke) :: tendency
546 real,
dimension(SZI_(G),SZJ_(G)) :: tendency_2d
547 real,
dimension(SZIB_(G),SZJ_(G)) :: trans_x_2d
548 real,
dimension(SZI_(G),SZJB_(G)) :: trans_y_2d
549 real,
dimension(G%ke) :: dtracer
553 integer :: i, j, k, m, ks, nk
555 real :: h_neglect, h_neglect_edge
557 if (.not.cs%remap_answers_2018)
then
558 h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
560 h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
570 if (tracer%id_dfxy_conc > 0 .or. tracer%id_dfxy_cont > 0 .or. tracer%id_dfxy_cont_2d > 0 .or. &
571 tracer%id_dfx_2d > 0 .or. tracer%id_dfy_2d > 0)
then
573 tendency(:,:,:) = 0.0
580 do j = g%jsc,g%jec ;
do i = g%isc-1,g%iec
581 if (g%mask2dCu(i,j)>0.)
then
582 call neutral_surface_flux(nk, cs%nsurf, cs%deg, h(i,j,:), h(i+1,j,:), &
583 tracer%t(i,j,:), tracer%t(i+1,j,:), &
584 cs%uPoL(i,j,:), cs%uPoR(i,j,:), &
585 cs%uKoL(i,j,:), cs%uKoR(i,j,:), &
586 cs%uhEff(i,j,:), uflx(i,j,:), &
587 cs%continuous_reconstruction, h_neglect, cs%remap_CS, h_neglect_edge)
592 do j = g%jsc-1,g%jec ;
do i = g%isc,g%iec
593 if (g%mask2dCv(i,j)>0.)
then
594 call neutral_surface_flux(nk, cs%nsurf, cs%deg, h(i,j,:), h(i,j+1,:), &
595 tracer%t(i,j,:), tracer%t(i,j+1,:), &
596 cs%vPoL(i,j,:), cs%vPoR(i,j,:), &
597 cs%vKoL(i,j,:), cs%vKoR(i,j,:), &
598 cs%vhEff(i,j,:), vflx(i,j,:), &
599 cs%continuous_reconstruction, h_neglect, cs%remap_CS, h_neglect_edge)
604 do j = g%jsc,g%jec ;
do i = g%isc,g%iec
605 if (g%mask2dT(i,j)>0.)
then
610 dtracer(k) = dtracer(k) + coef_x(i,j) * uflx(i,j,ks)
611 k = cs%uKoR(i-1,j,ks)
612 dtracer(k) = dtracer(k) - coef_x(i-1,j) * uflx(i-1,j,ks)
614 dtracer(k) = dtracer(k) + coef_y(i,j) * vflx(i,j,ks)
615 k = cs%vKoR(i,j-1,ks)
616 dtracer(k) = dtracer(k) - coef_y(i,j-1) * vflx(i,j-1,ks)
619 tracer%t(i,j,k) = tracer%t(i,j,k) + dtracer(k) * &
620 ( g%IareaT(i,j) / ( h(i,j,k) + gv%H_subroundoff ) )
623 if (tracer%id_dfxy_conc > 0 .or. tracer%id_dfxy_cont > 0 .or. tracer%id_dfxy_cont_2d > 0 )
then
625 tendency(i,j,k) = dtracer(k) * g%IareaT(i,j) * idt
634 if (tracer%id_dfx_2d > 0)
then
635 do j = g%jsc,g%jec ;
do i = g%isc-1,g%iec
637 if (g%mask2dCu(i,j)>0.)
then
639 trans_x_2d(i,j) = trans_x_2d(i,j) - coef_x(i,j) * uflx(i,j,ks)
641 trans_x_2d(i,j) = trans_x_2d(i,j) * idt
644 call post_data(tracer%id_dfx_2d, trans_x_2d(:,:), cs%diag)
649 if (tracer%id_dfy_2d > 0)
then
650 do j = g%jsc-1,g%jec ;
do i = g%isc,g%iec
652 if (g%mask2dCv(i,j)>0.)
then
654 trans_y_2d(i,j) = trans_y_2d(i,j) - coef_y(i,j) * vflx(i,j,ks)
656 trans_y_2d(i,j) = trans_y_2d(i,j) * idt
659 call post_data(tracer%id_dfy_2d, trans_y_2d(:,:), cs%diag)
663 if (tracer%id_dfxy_cont > 0)
then
664 call post_data(tracer%id_dfxy_cont, tendency(:,:,:), cs%diag)
668 if (tracer%id_dfxy_cont_2d > 0)
then
669 tendency_2d(:,:) = 0.
670 do j = g%jsc,g%jec ;
do i = g%isc,g%iec
672 tendency_2d(i,j) = tendency_2d(i,j) + tendency(i,j,k)
675 call post_data(tracer%id_dfxy_cont_2d, tendency_2d(:,:), cs%diag)
681 if (tracer%id_dfxy_conc > 0)
then
682 do k = 1, gv%ke ;
do j = g%jsc,g%jec ;
do i = g%isc,g%iec
683 tendency(i,j,k) = tendency(i,j,k) / ( h(i,j,k) + gv%H_subroundoff )
684 enddo ;
enddo ;
enddo
685 call post_data(tracer%id_dfxy_conc, tendency, cs%diag)
689 end subroutine neutral_diffusion
692 subroutine interface_scalar(nk, h, S, Si, i_method, h_neglect)
693 integer,
intent(in) :: nk
694 real,
dimension(nk),
intent(in) :: h
695 real,
dimension(nk),
intent(in) :: S
696 real,
dimension(nk+1),
intent(inout) :: Si
697 integer,
intent(in) :: i_method
699 real,
intent(in) :: h_neglect
701 integer :: k, km2, kp1
702 real,
dimension(nk) :: diff
705 call plm_diff(nk, h, s, 2, 1, diff)
706 si(1) = s(1) - 0.5 * diff(1)
707 if (i_method==1)
then
711 sa = s(k-1) + 0.5 * diff(k-1)
712 sb = s(k) - 0.5 * diff(k)
713 si(k) = 0.5 * ( sa + sb )
715 elseif (i_method==2)
then
721 si(k) = ppm_edge(h(km2), h(k-1), h(k), h(kp1), s(k-1), s(k), diff(k-1), diff(k), h_neglect)
724 si(nk+1) = s(nk) + 0.5 * diff(nk)
726 end subroutine interface_scalar
730 real function ppm_edge(hkm1, hk, hkp1, hkp2, Ak, Akp1, Pk, Pkp1, h_neglect)
731 real,
intent(in) :: hkm1
732 real,
intent(in) :: hk
733 real,
intent(in) :: hkp1
734 real,
intent(in) :: hkp2
735 real,
intent(in) :: ak
736 real,
intent(in) :: akp1
737 real,
intent(in) :: pk
738 real,
intent(in) :: pkp1
739 real,
intent(in) :: h_neglect
742 real :: r_hk_hkp1, r_2hk_hkp1, r_hk_2hkp1, f1, f2, f3, f4
744 r_hk_hkp1 = hk + hkp1
745 if (r_hk_hkp1 <= 0.)
then
746 ppm_edge = 0.5 * ( ak + akp1 )
749 r_hk_hkp1 = 1. / r_hk_hkp1
751 ppm_edge = ak + ( hk * r_hk_hkp1 ) * ( akp1 - ak )
753 ppm_edge = akp1 + ( hkp1 * r_hk_hkp1 ) * ( ak - akp1 )
756 r_2hk_hkp1 = 1. / ( ( 2. * hk + hkp1 ) + h_neglect )
757 r_hk_2hkp1 = 1. / ( ( hk + 2. * hkp1 ) + h_neglect )
758 f1 = 1./ ( ( hk + hkp1) + ( hkm1 + hkp2 ) )
759 f2 = 2. * ( hkp1 * hk ) * r_hk_hkp1 * &
760 ( ( hkm1 + hk ) * r_2hk_hkp1 - ( hkp2 + hkp1 ) * r_hk_2hkp1 )
761 f3 = hk * ( hkm1 + hk ) * r_2hk_hkp1
762 f4 = hkp1 * ( hkp1 + hkp2 ) * r_hk_2hkp1
764 ppm_edge = ppm_edge + f1 * ( f2 * ( akp1 - ak ) - ( f3 * pkp1 - f4 * pk ) )
766 end function ppm_edge
770 real function ppm_ave(xL, xR, aL, aR, aMean)
771 real,
intent(in) :: xl
772 real,
intent(in) :: xr
773 real,
intent(in) :: al
774 real,
intent(in) :: ar
775 real,
intent(in) :: amean
778 real :: dx, xave, a6, a6o3
781 xave = 0.5 * ( xr + xl )
782 a6o3 = 2. * amean - ( al + ar )
786 stop
'ppm_ave: dx<0 should not happend!'
788 stop
'ppm_ave: dx>1 should not happend!'
790 ppm_ave = al + ( ar - al ) * xr + a6 * xr * ( 1. - xr )
792 ppm_ave = ( al + xave * ( ( ar - al ) + a6 ) ) - a6o3 * ( xr**2 + xr * xl + xl**2 )
797 real function signum(a,x)
798 real,
intent(in) :: a
799 real,
intent(in) :: x
802 if (x==0.) signum = 0.
808 subroutine plm_diff(nk, h, S, c_method, b_method, diff)
809 integer,
intent(in) :: nk
810 real,
dimension(nk),
intent(in) :: h
811 real,
dimension(nk),
intent(in) :: S
812 integer,
intent(in) :: c_method
813 integer,
intent(in) :: b_method
814 real,
dimension(nk),
intent(inout) :: diff
824 real :: hkm1, hk, hkp1, Skm1, Sk, Skp1, diff_l, diff_r, diff_c
831 if ( ( hkp1 + hk ) * ( hkm1 + hk ) > 0.)
then
835 if (c_method==1)
then
837 if ( hk + 0.5 * (hkm1 + hkp1) /= 0. )
then
838 diff_c = ( skp1 - skm1 ) * ( hk / ( hk + 0.5 * (hkm1 + hkp1) ) )
842 elseif (c_method==2)
then
844 diff_c = fv_diff(hkm1, hk, hkp1, skm1, sk, skp1)
845 elseif (c_method==3)
then
847 diff_c = hk * fvlsq_slope(hkm1, hk, hkp1, skm1, sk, skp1)
850 diff_l = 2. * ( sk - skm1 )
851 diff_r = 2. * ( skp1 - sk )
852 if ( signum(1., diff_l) * signum(1., diff_r) <= 0. )
then
855 diff(k) = sign( min( abs(diff_l), abs(diff_c), abs(diff_r) ), diff_c )
861 if (b_method==1)
then
864 elseif (b_method==2)
then
865 diff(1) = ( s(2) - s(1) ) * 2. * ( h(1) / ( h(1) + h(2) ) )
866 diff(nk) = s(nk) - s(nk-1) * 2. * ( h(nk) / ( h(nk-1) + h(nk) ) )
869 end subroutine plm_diff
875 real function fv_diff(hkm1, hk, hkp1, Skm1, Sk, Skp1)
876 real,
intent(in) :: hkm1
877 real,
intent(in) :: hk
878 real,
intent(in) :: hkp1
879 real,
intent(in) :: skm1
880 real,
intent(in) :: sk
881 real,
intent(in) :: skp1
884 real :: h_sum, hp, hm
886 h_sum = ( hkm1 + hkp1 ) + hk
887 if (h_sum /= 0.) h_sum = 1./ h_sum
889 if (hm /= 0.) hm = 1./ hm
891 if (hp /= 0.) hp = 1./ hp
892 fv_diff = ( hk * h_sum ) * &
893 ( ( 2. * hkm1 + hk ) * hp * ( skp1 - sk ) &
894 + ( 2. * hkp1 + hk ) * hm * ( sk - skm1 ) )
901 real function fvlsq_slope(hkm1, hk, hkp1, Skm1, Sk, Skp1)
902 real,
intent(in) :: hkm1
903 real,
intent(in) :: hk
904 real,
intent(in) :: hkp1
905 real,
intent(in) :: skm1
906 real,
intent(in) :: sk
907 real,
intent(in) :: skp1
911 real :: h_sum, hx_sum, hxsq_sum, hxy_sum, hy_sum, det
913 xkm1 = -0.5 * ( hk + hkm1 )
914 xkp1 = 0.5 * ( hk + hkp1 )
915 h_sum = ( hkm1 + hkp1 ) + hk
916 hx_sum = hkm1*xkm1 + hkp1*xkp1
917 hxsq_sum = hkm1*(xkm1**2) + hkp1*(xkp1**2)
918 hxy_sum = hkm1*xkm1*skm1 + hkp1*xkp1*skp1
919 hy_sum = ( hkm1*skm1 + hkp1*skp1 ) + hk*sk
920 det = h_sum * hxsq_sum - hx_sum**2
923 fvlsq_slope = ( h_sum * hxy_sum - hx_sum*hy_sum ) / det
927 end function fvlsq_slope
931 subroutine find_neutral_surface_positions_continuous(nk, Pl, Tl, Sl, dRdTl, dRdSl, Pr, Tr, Sr, &
932 dRdTr, dRdSr, PoL, PoR, KoL, KoR, hEff, bl_kl, bl_kr, bl_zl, bl_zr)
933 integer,
intent(in) :: nk
934 real,
dimension(nk+1),
intent(in) :: Pl
935 real,
dimension(nk+1),
intent(in) :: Tl
936 real,
dimension(nk+1),
intent(in) :: Sl
937 real,
dimension(nk+1),
intent(in) :: dRdTl
938 real,
dimension(nk+1),
intent(in) :: dRdSl
939 real,
dimension(nk+1),
intent(in) :: Pr
940 real,
dimension(nk+1),
intent(in) :: Tr
941 real,
dimension(nk+1),
intent(in) :: Sr
942 real,
dimension(nk+1),
intent(in) :: dRdTr
943 real,
dimension(nk+1),
intent(in) :: dRdSr
944 real,
dimension(2*nk+2),
intent(inout) :: PoL
946 real,
dimension(2*nk+2),
intent(inout) :: PoR
948 integer,
dimension(2*nk+2),
intent(inout) :: KoL
949 integer,
dimension(2*nk+2),
intent(inout) :: KoR
950 real,
dimension(2*nk+1),
intent(inout) :: hEff
952 integer,
optional,
intent(in) :: bl_kl
953 integer,
optional,
intent(in) :: bl_kr
954 real,
optional,
intent(in) :: bl_zl
955 real,
optional,
intent(in) :: bl_zr
963 logical :: searching_left_column
964 logical :: searching_right_column
965 logical :: reached_bottom
966 integer :: krm1, klm1
967 real :: dRho, dRhoTop, dRhoBot
969 integer :: lastK_left, lastK_right
970 real :: lastP_left, lastP_right
971 logical :: interior_limit
982 reached_bottom = .false.
985 interior_limit =
PRESENT(bl_kl) .and.
PRESENT(bl_kr) .and.
PRESENT(bl_zr) .and.
PRESENT(bl_zl)
988 neutral_surfaces:
do k_surface = 1, ns
990 if (klm1>nk) stop
'find_neutral_surface_positions(): klm1 went out of bounds!'
992 if (krm1>nk) stop
'find_neutral_surface_positions(): krm1 went out of bounds!'
995 drho = 0.5 * ( ( drdtr(kr) + drdtl(kl) ) * ( tr(kr) - tl(kl) ) &
996 + ( drdsr(kr) + drdsl(kl) ) * ( sr(kr) - sl(kl) ) )
998 if (.not. reached_bottom)
then
1000 searching_left_column = .true.
1001 searching_right_column = .false.
1002 elseif (drho > 0.)
then
1003 searching_right_column = .true.
1004 searching_left_column = .false.
1006 if (kl + kr == 2)
then
1007 searching_left_column = .true.
1008 searching_right_column = .false.
1010 searching_left_column = .not. searching_left_column
1011 searching_right_column = .not. searching_right_column
1016 if (searching_left_column)
then
1019 drhotop = 0.5 * ( ( drdtl(klm1) + drdtr(kr) ) * ( tl(klm1) - tr(kr) ) &
1020 + ( drdsl(klm1) + drdsr(kr) ) * ( sl(klm1) - sr(kr) ) )
1022 drhobot = 0.5 * ( ( drdtl(klm1+1) + drdtr(kr) ) * ( tl(klm1+1) - tr(kr) ) &
1023 + ( drdsl(klm1+1) + drdsr(kr) ) * ( sl(klm1+1) - sr(kr) ) )
1027 if (drhotop > 0. .or. kr+kl==2)
then
1029 elseif (drhotop >= drhobot)
then
1034 pol(k_surface) = interpolate_for_nondim_position( drhotop, pl(klm1), drhobot, pl(klm1+1) )
1036 if (pol(k_surface)>=1. .and. klm1<nk)
then
1038 pol(k_surface) = pol(k_surface) - 1.
1040 if (real(klm1-lastk_left)+(pol(k_surface)-lastp_left)<0.)
then
1041 pol(k_surface) = lastp_left
1044 kol(k_surface) = klm1
1055 reached_bottom = .true.
1056 searching_right_column = .true.
1057 searching_left_column = .false.
1059 elseif (searching_right_column)
then
1062 drhotop = 0.5 * ( ( drdtr(krm1) + drdtl(kl) ) * ( tr(krm1) - tl(kl) ) + &
1063 ( drdsr(krm1) + drdsl(kl) ) * ( sr(krm1) - sl(kl) ) )
1065 drhobot = 0.5 * ( ( drdtr(krm1+1) + drdtl(kl) ) * ( tr(krm1+1) - tl(kl) ) + &
1066 ( drdsr(krm1+1) + drdsl(kl) ) * ( sr(krm1+1) - sl(kl) ) )
1070 if (drhotop >= 0. .or. kr+kl==2)
then
1072 elseif (drhotop >= drhobot)
then
1077 por(k_surface) = interpolate_for_nondim_position( drhotop, pr(krm1), drhobot, pr(krm1+1) )
1079 if (por(k_surface)>=1. .and. krm1<nk)
then
1081 por(k_surface) = por(k_surface) - 1.
1083 if (real(krm1-lastk_right)+(por(k_surface)-lastp_right)<0.)
then
1084 por(k_surface) = lastp_right
1087 kor(k_surface) = krm1
1098 reached_bottom = .true.
1099 searching_right_column = .false.
1100 searching_left_column = .true.
1105 if (interior_limit)
then
1106 if (kol(k_surface)<=bl_kl)
then
1107 kol(k_surface) = bl_kl
1108 if (pol(k_surface)<bl_zl)
then
1109 pol(k_surface) = bl_zl
1112 if (kor(k_surface)<=bl_kr)
then
1113 kor(k_surface) = bl_kr
1114 if (por(k_surface)<bl_zr)
then
1115 por(k_surface) = bl_zr
1120 lastk_left = kol(k_surface) ; lastp_left = pol(k_surface)
1121 lastk_right = kor(k_surface) ; lastp_right = por(k_surface)
1125 if (k_surface>1)
then
1126 hl = absolute_position(nk,ns,pl,kol,pol,k_surface) - absolute_position(nk,ns,pl,kol,pol,k_surface-1)
1127 hr = absolute_position(nk,ns,pr,kor,por,k_surface) - absolute_position(nk,ns,pr,kor,por,k_surface-1)
1128 if ( hl + hr > 0.)
then
1129 heff(k_surface-1) = 2. * hl * hr / ( hl + hr )
1131 heff(k_surface-1) = 0.
1135 enddo neutral_surfaces
1137 end subroutine find_neutral_surface_positions_continuous
1142 real function interpolate_for_nondim_position(dRhoNeg, Pneg, dRhoPos, Ppos)
1143 real,
intent(in) :: drhoneg
1144 real,
intent(in) :: pneg
1145 real,
intent(in) :: drhopos
1146 real,
intent(in) :: ppos
1148 character(len=120) :: mesg
1150 if (ppos < pneg)
then
1151 call mom_error(fatal,
'interpolate_for_nondim_position: Houston, we have a problem! Ppos<Pneg')
1152 elseif (drhoneg>drhopos)
then
1153 write(stderr,*)
'dRhoNeg, Pneg, dRhoPos, Ppos=',drhoneg, pneg, drhopos, ppos
1154 write(mesg,*)
'dRhoNeg, Pneg, dRhoPos, Ppos=', drhoneg, pneg, drhopos, ppos
1155 call mom_error(warning,
'interpolate_for_nondim_position: '//trim(mesg))
1156 elseif (drhoneg>drhopos)
then
1157 call mom_error(fatal,
'interpolate_for_nondim_position: Houston, we have a problem! dRhoNeg>dRhoPos')
1159 if (ppos<=pneg)
then
1160 interpolate_for_nondim_position = 0.5
1161 elseif ( drhopos - drhoneg > 0. )
then
1162 interpolate_for_nondim_position = min( 1., max( 0., -drhoneg / ( drhopos - drhoneg ) ) )
1163 elseif ( drhopos - drhoneg == 0)
then
1164 if (drhoneg>0.)
then
1165 interpolate_for_nondim_position = 0.
1166 elseif (drhoneg<0.)
then
1167 interpolate_for_nondim_position = 1.
1169 interpolate_for_nondim_position = 0.5
1172 interpolate_for_nondim_position = 0.5
1174 if ( interpolate_for_nondim_position < 0. ) &
1175 call mom_error(fatal,
'interpolate_for_nondim_position: Houston, we have a problem! Pint < Pneg')
1176 if ( interpolate_for_nondim_position > 1. ) &
1177 call mom_error(fatal,
'interpolate_for_nondim_position: Houston, we have a problem! Pint > Ppos')
1178 end function interpolate_for_nondim_position
1183 subroutine find_neutral_surface_positions_discontinuous(CS, nk, &
1184 Pres_l, hcol_l, Tl, Sl, ppoly_T_l, ppoly_S_l, stable_l, &
1185 Pres_r, hcol_r, Tr, Sr, ppoly_T_r, ppoly_S_r, stable_r, &
1186 PoL, PoR, KoL, KoR, hEff, zeta_bot_L, zeta_bot_R, k_bot_L, k_bot_R, hard_fail_heff)
1189 integer,
intent(in) :: nk
1190 real,
dimension(nk,2),
intent(in) :: Pres_l
1191 real,
dimension(nk),
intent(in) :: hcol_l
1193 real,
dimension(nk,2),
intent(in) :: Tl
1194 real,
dimension(nk,2),
intent(in) :: Sl
1195 real,
dimension(:,:),
intent(in) :: ppoly_T_l
1196 real,
dimension(:,:),
intent(in) :: ppoly_S_l
1197 logical,
dimension(nk),
intent(in) :: stable_l
1198 real,
dimension(nk,2),
intent(in) :: Pres_r
1199 real,
dimension(nk),
intent(in) :: hcol_r
1201 real,
dimension(nk,2),
intent(in) :: Tr
1202 real,
dimension(nk,2),
intent(in) :: Sr
1203 real,
dimension(:,:),
intent(in) :: ppoly_T_r
1204 real,
dimension(:,:),
intent(in) :: ppoly_S_r
1205 logical,
dimension(nk),
intent(in) :: stable_r
1206 real,
dimension(4*nk),
intent(inout) :: PoL
1208 real,
dimension(4*nk),
intent(inout) :: PoR
1210 integer,
dimension(4*nk),
intent(inout) :: KoL
1211 integer,
dimension(4*nk),
intent(inout) :: KoR
1212 real,
dimension(4*nk-1),
intent(inout) :: hEff
1214 real,
optional,
intent(in) :: zeta_bot_L
1216 real,
optional,
intent(in) :: zeta_bot_R
1219 integer,
optional,
intent(in) :: k_bot_L
1220 integer,
optional,
intent(in) :: k_bot_R
1221 logical,
optional,
intent(in) :: hard_fail_heff
1225 integer :: k_surface
1226 integer :: kl_left, kl_right
1227 integer :: ki_left, ki_right
1228 logical :: searching_left_column
1229 logical :: searching_right_column
1230 logical :: reached_bottom
1231 logical :: search_layer
1232 logical :: fail_heff
1236 real :: lastP_left, lastP_right
1237 integer :: k_init_L, k_init_R
1238 real :: p_init_L, p_init_R
1247 reached_bottom = .false.
1248 searching_left_column = .false.
1249 searching_right_column = .false.
1252 if (
PRESENT(hard_fail_heff)) fail_heff = hard_fail_heff
1254 if (
PRESENT(k_bot_l) .and.
PRESENT(k_bot_r) .and.
PRESENT(zeta_bot_l) .and.
PRESENT(zeta_bot_r))
then
1255 k_init_l = k_bot_l; k_init_r = k_bot_r
1256 p_init_l = zeta_bot_l; p_init_r = zeta_bot_r
1257 lastp_left = zeta_bot_l; lastp_right = zeta_bot_r
1258 kl_left = k_bot_l; kl_right = k_bot_r
1260 k_init_l = 1 ; k_init_r = 1
1261 p_init_l = 0. ; p_init_r = 0.
1264 neutral_surfaces:
do k_surface = 1, ns
1266 if (k_surface == ns)
then
1272 elseif (.not. stable_l(kl_left))
then
1273 if (k_surface > 1)
then
1274 pol(k_surface) = ki_left - 1
1275 kol(k_surface) = kl_left
1276 por(k_surface) = por(k_surface-1)
1277 kor(k_surface) = kor(k_surface-1)
1279 por(k_surface) = p_init_r
1280 kor(k_surface) = k_init_r
1281 pol(k_surface) = p_init_l
1282 kol(k_surface) = k_init_l
1284 call increment_interface(nk, kl_left, ki_left, reached_bottom, searching_left_column, searching_right_column)
1285 searching_left_column = .true.
1286 searching_right_column = .false.
1287 elseif (.not. stable_r(kl_right))
then
1288 if (k_surface > 1)
then
1289 por(k_surface) = ki_right - 1
1290 kor(k_surface) = kl_right
1291 pol(k_surface) = pol(k_surface-1)
1292 kol(k_surface) = kol(k_surface-1)
1299 call increment_interface(nk, kl_right, ki_right, reached_bottom, searching_right_column, searching_left_column)
1300 searching_left_column = .false.
1301 searching_right_column = .true.
1306 call calc_delta_rho_and_derivs(cs, &
1307 tr(kl_right, ki_right), sr(kl_right, ki_right), pres_r(kl_right,ki_right), &
1308 tl(kl_left, ki_left), sl(kl_left, ki_left) , pres_l(kl_left,ki_left), &
1310 if (cs%debug)
write(stdout,
'(A,I2,A,E12.4,A,I2,A,I2,A,I2,A,I2)') &
1311 "k_surface=",k_surface,
" dRho=",cs%R_to_kg_m3*drho, &
1312 "kl_left=",kl_left,
" ki_left=",ki_left,
" kl_right=",kl_right,
" ki_right=",ki_right
1314 if (.not. reached_bottom)
then
1316 searching_left_column = .true.
1317 searching_right_column = .false.
1318 elseif (drho > 0.)
then
1319 searching_left_column = .false.
1320 searching_right_column = .true.
1322 if ( ( kl_left + kl_right == 2 ) .and. (ki_left + ki_right == 2) )
then
1323 searching_left_column = .true.
1324 searching_right_column = .false.
1326 searching_left_column = .not. searching_left_column
1327 searching_right_column = .not. searching_right_column
1331 if (searching_left_column)
then
1333 por(k_surface) = ki_right - 1.
1334 kor(k_surface) = kl_right
1335 pol(k_surface) = search_other_column(cs, k_surface, lastp_left, &
1336 tr(kl_right, ki_right), sr(kl_right, ki_right), pres_r(kl_right, ki_right), &
1337 tl(kl_left,1), sl(kl_left,1), pres_l(kl_left,1), &
1338 tl(kl_left,2), sl(kl_left,2), pres_l(kl_left,2), &
1339 ppoly_t_l(kl_left,:), ppoly_s_l(kl_left,:))
1340 kol(k_surface) = kl_left
1343 write(stdout,
'(A,I2)')
"Searching left layer ", kl_left
1344 write(stdout,
'(A,I2,X,I2)')
"Searching from right: ", kl_right, ki_right
1345 write(stdout,*)
"Temp/Salt Reference: ", tr(kl_right,ki_right), sr(kl_right,ki_right)
1346 write(stdout,*)
"Temp/Salt Top L: ", tl(kl_left,1), sl(kl_left,1)
1347 write(stdout,*)
"Temp/Salt Bot L: ", tl(kl_left,2), sl(kl_left,2)
1349 call increment_interface(nk, kl_right, ki_right, reached_bottom, searching_right_column, searching_left_column)
1350 lastp_left = pol(k_surface)
1352 if ( kl_right == (kor(k_surface) + 1) ) lastp_right = 0.
1354 elseif (searching_right_column)
then
1356 pol(k_surface) = ki_left - 1.
1357 kol(k_surface) = kl_left
1358 por(k_surface) = search_other_column(cs, k_surface, lastp_right, &
1359 tl(kl_left, ki_left), sl(kl_left, ki_left), pres_l(kl_left, ki_left), &
1360 tr(kl_right,1), sr(kl_right,1), pres_r(kl_right,1), &
1361 tr(kl_right,2), sr(kl_right,2), pres_r(kl_right,2), &
1362 ppoly_t_r(kl_right,:), ppoly_s_r(kl_right,:))
1363 kor(k_surface) = kl_right
1366 write(stdout,
'(A,I2)')
"Searching right layer ", kl_right
1367 write(stdout,
'(A,I2,X,I2)')
"Searching from left: ", kl_left, ki_left
1368 write(stdout,*)
"Temp/Salt Reference: ", tl(kl_left,ki_left), sl(kl_left,ki_left)
1369 write(stdout,*)
"Temp/Salt Top L: ", tr(kl_right,1), sr(kl_right,1)
1370 write(stdout,*)
"Temp/Salt Bot L: ", tr(kl_right,2), sr(kl_right,2)
1372 call increment_interface(nk, kl_left, ki_left, reached_bottom, searching_left_column, searching_right_column)
1373 lastp_right = por(k_surface)
1375 if ( kl_left == (kol(k_surface) + 1) ) lastp_left = 0.
1379 if (cs%debug)
write(stdout,
'(A,I3,A,ES16.6,A,I2,A,ES16.6)')
"KoL:", kol(k_surface),
" PoL:", pol(k_surface), &
1380 " KoR:", kor(k_surface),
" PoR:", por(k_surface)
1383 if (k_surface>1)
then
1384 if ( kol(k_surface) == kol(k_surface-1) .and. kor(k_surface) == kor(k_surface-1) )
then
1385 hl = (pol(k_surface) - pol(k_surface-1))*hcol_l(kol(k_surface))
1386 hr = (por(k_surface) - por(k_surface-1))*hcol_r(kor(k_surface))
1387 if (hl < 0. .or. hr < 0.)
then
1389 call mom_error(fatal,
"Negative thicknesses in neutral diffusion")
1391 if (searching_left_column)
then
1392 pol(k_surface) = pol(k_surface-1)
1393 kol(k_surface) = kol(k_surface-1)
1394 elseif (searching_right_column)
then
1395 por(k_surface) = por(k_surface-1)
1396 kor(k_surface) = kor(k_surface-1)
1399 elseif ( hl + hr == 0. )
then
1400 heff(k_surface-1) = 0.
1402 heff(k_surface-1) = 2. * ( (hl * hr) / ( hl + hr ) )
1403 if ( kol(k_surface) /= kol(k_surface-1) )
then
1404 call mom_error(fatal,
"Neutral sublayer spans multiple layers")
1406 if ( kor(k_surface) /= kor(k_surface-1) )
then
1407 call mom_error(fatal,
"Neutral sublayer spans multiple layers")
1411 heff(k_surface-1) = 0.
1414 enddo neutral_surfaces
1415 end subroutine find_neutral_surface_positions_discontinuous
1418 subroutine mark_unstable_cells(CS, nk, T, S, P, stable_cell)
1420 integer,
intent(in) :: nk
1421 real,
dimension(nk,2),
intent(in) :: T
1422 real,
dimension(nk,2),
intent(in) :: S
1423 real,
dimension(nk,2),
intent(in) :: P
1424 logical,
dimension(nk),
intent( out) :: stable_cell
1426 integer :: k, first_stable, prev_stable
1430 call calc_delta_rho_and_derivs( cs, t(k,2), s(k,2), max(p(k,2), cs%ref_pres), &
1431 t(k,1), s(k,1), max(p(k,1), cs%ref_pres), delta_rho )
1432 stable_cell(k) = (delta_rho > 0.)
1434 end subroutine mark_unstable_cells
1437 real function search_other_column(CS, ksurf, pos_last, T_from, S_from, P_from, T_top, S_top, P_top, &
1438 T_bot, S_bot, P_bot, T_poly, S_poly )
result(pos)
1440 integer,
intent(in ) :: ksurf
1441 real,
intent(in ) :: pos_last
1443 real,
intent(in ) :: t_from
1444 real,
intent(in ) :: s_from
1445 real,
intent(in ) :: p_from
1446 real,
intent(in ) :: t_top
1447 real,
intent(in ) :: s_top
1448 real,
intent(in ) :: p_top
1450 real,
intent(in ) :: t_bot
1451 real,
intent(in ) :: s_bot
1452 real,
intent(in ) :: p_bot
1454 real,
dimension(:),
intent(in ) :: t_poly
1455 real,
dimension(:),
intent(in ) :: s_poly
1457 real :: drhotop, drhobot
1458 real :: drdt_top, drdt_bot, drdt_from
1459 real :: drds_top, drds_bot, drds_from
1462 if (cs%neutral_pos_method == 1 .or. cs%neutral_pos_method == 3)
then
1463 call calc_delta_rho_and_derivs(cs, t_top, s_top, p_top, t_from, s_from, p_from, drhotop)
1464 call calc_delta_rho_and_derivs(cs, t_bot, s_bot, p_bot, t_from, s_from, p_from, drhobot)
1465 elseif (cs%neutral_pos_method == 2)
then
1466 call calc_delta_rho_and_derivs(cs, t_top, s_top, p_top, t_from, s_from, p_from, drhotop, &
1467 drdt_top, drds_top, drdt_from, drds_from)
1468 call calc_delta_rho_and_derivs(cs, t_bot, s_bot, p_bot, t_from, s_from, p_from, drhobot, &
1469 drdt_bot, drds_bot, drdt_from, drds_from)
1473 if ( (drhotop > 0.) .or. (ksurf == 1) )
then
1475 elseif ( drhotop > drhobot )
then
1477 elseif ( drhotop < 0. .and. drhobot < 0.)
then
1479 elseif ( drhotop == 0. .and. drhobot == 0. )
then
1481 elseif ( drhobot == 0. )
then
1483 elseif ( drhotop == 0. )
then
1492 if (cs%neutral_pos_method==1)
then
1493 pos = interpolate_for_nondim_position( drhotop, p_top, drhobot, p_bot )
1496 elseif (cs%neutral_pos_method == 2)
then
1497 pos = find_neutral_pos_linear( cs, pos_last, t_from, s_from, drdt_from, drds_from, &
1498 drdt_top, drds_top, drdt_bot, drds_bot, t_poly, s_poly )
1499 elseif (cs%neutral_pos_method == 3)
then
1500 pos = find_neutral_pos_full( cs, pos_last, t_from, s_from, p_from, p_top, p_bot, t_poly, s_poly)
1503 end function search_other_column
1506 subroutine increment_interface(nk, kl, ki, reached_bottom, searching_this_column, searching_other_column)
1507 integer,
intent(in ) :: nk
1508 integer,
intent(inout) :: kl
1509 integer,
intent(inout) :: ki
1510 logical,
intent(inout) :: reached_bottom
1511 logical,
intent(inout) :: searching_this_column
1512 logical,
intent(inout) :: searching_other_column
1515 reached_bottom = .false.
1517 if ((ki == 2) .and. (kl < nk) )
then
1520 elseif ((kl == nk) .and. (ki==2))
then
1521 reached_bottom = .true.
1522 searching_this_column = .false.
1523 searching_other_column = .true.
1528 call mom_error(fatal,
"Unanticipated eventuality in increment_interface")
1530 end subroutine increment_interface
1540 function find_neutral_pos_linear( CS, z0, T_ref, S_ref, dRdT_ref, dRdS_ref, &
1541 dRdT_top, dRdS_top, dRdT_bot, dRdS_bot, ppoly_T, ppoly_S )
result( z )
1543 real,
intent(in) :: z0
1545 real,
intent(in) :: t_ref
1546 real,
intent(in) :: s_ref
1547 real,
intent(in) :: drdt_ref
1549 real,
intent(in) :: drds_ref
1551 real,
intent(in) :: drdt_top
1553 real,
intent(in) :: drds_top
1555 real,
intent(in) :: drdt_bot
1557 real,
intent(in) :: drds_bot
1559 real,
dimension(:),
intent(in) :: ppoly_t
1561 real,
dimension(:),
intent(in) :: ppoly_s
1569 real :: drho, drho_dz
1574 real :: drho_min, drho_max
1575 real :: ztest, zmin, zmax
1581 nterm =
SIZE(ppoly_t)
1584 drdt_diff = drdt_bot - drdt_top
1585 drds_diff = drds_bot - drds_top
1591 t_z = evaluation_polynomial( ppoly_t, nterm, zmin )
1592 s_z = evaluation_polynomial( ppoly_s, nterm, zmin )
1593 drdt_z = a1*drdt_top + a2*drdt_bot
1594 drds_z = a1*drds_top + a2*drds_bot
1595 drho_min = 0.5*((drdt_z+drdt_ref)*(t_z-t_ref) + (drds_z+drds_ref)*(s_z-s_ref))
1597 t_z = evaluation_polynomial( ppoly_t, nterm, 1. )
1598 s_z = evaluation_polynomial( ppoly_s, nterm, 1. )
1599 drho_max = 0.5*((drdt_bot+drdt_ref)*(t_z-t_ref) + (drds_bot+drds_ref)*(s_z-s_ref))
1601 if (drho_min >= 0.)
then
1604 elseif (drho_max == 0.)
then
1608 if ( sign(1.,drho_min) == sign(1.,drho_max) )
then
1609 call mom_error(fatal,
"drho_min is the same sign as dhro_max")
1614 do iter = 1, cs%max_iter
1618 drdt_z = a1*drdt_top + a2*drdt_bot
1619 drds_z = a1*drds_top + a2*drds_bot
1620 t_z = evaluation_polynomial( ppoly_t, nterm, z )
1621 s_z = evaluation_polynomial( ppoly_s, nterm, z )
1622 drho = 0.5*((drdt_z+drdt_ref)*(t_z-t_ref) + (drds_z+drds_ref)*(s_z-s_ref))
1625 if (abs(drho) <= cs%drho_tol)
exit
1627 if (drho < 0. .and. drho > drho_min)
then
1630 elseif (drho > 0. .and. drho < drho_max)
then
1636 dt_dz = first_derivative_polynomial( ppoly_t, nterm, z )
1637 ds_dz = first_derivative_polynomial( ppoly_s, nterm, z )
1638 drho_dz = 0.5*( (drdt_diff*(t_z - t_ref) + (drdt_ref+drdt_z)*dt_dz) + &
1639 (drds_diff*(s_z - s_ref) + (drds_ref+drds_z)*ds_dz) )
1641 ztest = z - drho/drho_dz
1643 if (ztest < zmin .or. ztest > zmax)
then
1644 if ( drho < 0. )
then
1645 ztest = 0.5*(z + zmax)
1647 ztest = 0.5*(zmin + z)
1652 if ( abs(z-ztest) <= cs%x_tol )
exit
1657 end function find_neutral_pos_linear
1661 function find_neutral_pos_full( CS, z0, T_ref, S_ref, P_ref, P_top, P_bot, ppoly_T, ppoly_S )
result( z )
1663 real,
intent(in) :: z0
1665 real,
intent(in) :: t_ref
1666 real,
intent(in) :: s_ref
1667 real,
intent(in) :: p_ref
1668 real,
intent(in) :: p_top
1669 real,
intent(in) :: p_bot
1670 real,
dimension(:),
intent(in) :: ppoly_t
1672 real,
dimension(:),
intent(in) :: ppoly_s
1679 real :: drho_a, drho_b, drho_c
1689 nterm =
SIZE(ppoly_t)
1692 tb = evaluation_polynomial( ppoly_t, nterm, b )
1693 sb = evaluation_polynomial( ppoly_s, nterm, b )
1694 pb = p_top*(1.-b) + p_bot*b
1695 call calc_delta_rho_and_derivs(cs, tb, sb, pb, t_ref, s_ref, p_ref, drho_b)
1698 tc = evaluation_polynomial( ppoly_t, nterm, 1. )
1699 sc = evaluation_polynomial( ppoly_s, nterm, 1. )
1701 call calc_delta_rho_and_derivs(cs, tc, sc, pc, t_ref, s_ref, p_ref, drho_c)
1703 if (drho_b >= 0.)
then
1706 elseif (drho_c == 0.)
then
1710 if ( sign(1.,drho_b) == sign(1.,drho_c) )
then
1715 do iter = 1, cs%max_iter
1717 a = (drho_b*c - drho_c*b)/(drho_b-drho_c)
1718 ta = evaluation_polynomial( ppoly_t, nterm, a )
1719 sa = evaluation_polynomial( ppoly_s, nterm, a )
1720 pa = p_top*(1.-a) + p_bot*a
1721 call calc_delta_rho_and_derivs(cs, ta, sa, pa, t_ref, s_ref, p_ref, drho_a)
1722 if (abs(drho_a) < cs%drho_tol)
then
1727 if (drho_a*drho_c > 0.)
then
1728 if ( abs(a-c)<cs%x_tol)
then
1732 c = a ; drho_c = drho_a;
1733 if (side == -1) drho_b = 0.5*drho_b
1735 elseif ( drho_b*drho_a > 0 )
then
1736 if ( abs(a-b)<cs%x_tol)
then
1740 b = a ; drho_b = drho_a
1741 if (side == 1) drho_c = 0.5*drho_c
1751 end function find_neutral_pos_full
1754 subroutine calc_delta_rho_and_derivs(CS, T1, S1, p1_in, T2, S2, p2_in, drho, &
1755 drdt1_out, drds1_out, drdt2_out, drds2_out )
1757 real,
intent(in ) :: T1
1758 real,
intent(in ) :: S1
1759 real,
intent(in ) :: p1_in
1760 real,
intent(in ) :: T2
1761 real,
intent(in ) :: S2
1762 real,
intent(in ) :: p2_in
1763 real,
intent( out) :: drho
1764 real,
optional,
intent( out) :: dRdT1_out
1765 real,
optional,
intent( out) :: dRdS1_out
1766 real,
optional,
intent( out) :: dRdT2_out
1767 real,
optional,
intent( out) :: dRdS2_out
1770 real :: p1, p2, pmid
1771 real :: drdt1, drdt2
1772 real :: drds1, drds2
1773 real :: drdp1, drdp2
1776 if (cs%ref_pres > 0.)
then
1785 if (trim(cs%delta_rho_form) ==
'full')
then
1786 pmid = 0.5 * (p1 + p2)
1791 elseif (trim(cs%delta_rho_form) ==
'mid_pressure')
then
1792 pmid = 0.5 * (p1 + p2)
1793 if (cs%ref_pres>=0) pmid = cs%ref_pres
1796 drho = delta_rho_from_derivs( t1, s1, p1, drdt1, drds1, t2, s2, p2, drdt2, drds2)
1797 elseif (trim(cs%delta_rho_form) ==
'local_pressure')
then
1800 drho = delta_rho_from_derivs( t1, s1, p1, drdt1, drds1, t2, s2, p2, drdt2, drds2)
1802 call mom_error(fatal,
"delta_rho_form is not recognized")
1805 if (
PRESENT(drdt1_out)) drdt1_out = drdt1
1806 if (
PRESENT(drds1_out)) drds1_out = drds1
1807 if (
PRESENT(drdt2_out)) drdt2_out = drdt2
1808 if (
PRESENT(drds2_out)) drds2_out = drds2
1810 end subroutine calc_delta_rho_and_derivs
1816 function delta_rho_from_derivs( T1, S1, P1, dRdT1, dRdS1, &
1817 T2, S2, P2, dRdT2, dRdS2 )
result (drho)
1831 drho = 0.5 * ( (drdt1+drdt2)*(t1-t2) + (drds1+drds2)*(s1-s2))
1833 end function delta_rho_from_derivs
1836 function absolute_position(n,ns,Pint,Karr,NParr,k_surface)
1837 integer,
intent(in) :: n
1838 integer,
intent(in) :: ns
1839 real,
intent(in) :: pint(n+1)
1840 integer,
intent(in) :: karr(ns)
1841 real,
intent(in) :: nparr(ns)
1842 integer,
intent(in) :: k_surface
1843 real :: absolute_position
1849 if (k>n) stop
'absolute_position: k>nk is out of bounds!'
1850 absolute_position = pint(k) + nparr(k_surface) * ( pint(k+1) - pint(k) )
1852 end function absolute_position
1855 function absolute_positions(n,ns,Pint,Karr,NParr)
1856 integer,
intent(in) :: n
1857 integer,
intent(in) :: ns
1858 real,
intent(in) :: pint(n+1)
1859 integer,
intent(in) :: karr(ns)
1860 real,
intent(in) :: nparr(ns)
1862 real,
dimension(ns) :: absolute_positions
1866 integer :: k_surface, k
1868 do k_surface = 1, ns
1869 absolute_positions(k_surface) = absolute_position(n,ns,pint,karr,nparr,k_surface)
1872 end function absolute_positions
1875 subroutine neutral_surface_flux(nk, nsurf, deg, hl, hr, Tl, Tr, PiL, PiR, KoL, KoR, &
1876 hEff, Flx, continuous, h_neglect, remap_CS, h_neglect_edge)
1877 integer,
intent(in) :: nk
1878 integer,
intent(in) :: nsurf
1879 integer,
intent(in) :: deg
1880 real,
dimension(nk),
intent(in) :: hl
1881 real,
dimension(nk),
intent(in) :: hr
1882 real,
dimension(nk),
intent(in) :: Tl
1883 real,
dimension(nk),
intent(in) :: Tr
1884 real,
dimension(nsurf),
intent(in) :: PiL
1886 real,
dimension(nsurf),
intent(in) :: PiR
1888 integer,
dimension(nsurf),
intent(in) :: KoL
1889 integer,
dimension(nsurf),
intent(in) :: KoR
1890 real,
dimension(nsurf-1),
intent(in) :: hEff
1892 real,
dimension(nsurf-1),
intent(inout) :: Flx
1893 logical,
intent(in) :: continuous
1894 real,
intent(in) :: h_neglect
1898 real,
optional,
intent(in) :: h_neglect_edge
1901 integer :: k_sublayer, klb, klt, krb, krt, k
1902 real :: T_right_top, T_right_bottom, T_right_layer, T_right_sub, T_right_top_int, T_right_bot_int
1903 real :: T_left_top, T_left_bottom, T_left_layer, T_left_sub, T_left_top_int, T_left_bot_int
1904 real :: dT_top, dT_bottom, dT_layer, dT_ave, dT_sublayer, dT_top_int, dT_bot_int
1905 real,
dimension(nk+1) :: Til
1906 real,
dimension(nk+1) :: Tir
1907 real,
dimension(nk) :: aL_l
1908 real,
dimension(nk) :: aR_l
1909 real,
dimension(nk) :: aL_r
1910 real,
dimension(nk) :: aR_r
1913 real,
dimension(nk,2) :: Tid_l
1914 real,
dimension(nk,2) :: Tid_r
1915 real,
dimension(nk,deg+1) :: ppoly_r_coeffs_l
1916 real,
dimension(nk,deg+1) :: ppoly_r_coeffs_r
1917 real,
dimension(nk,deg+1) :: ppoly_r_S_l
1918 real,
dimension(nk,deg+1) :: ppoly_r_S_r
1919 logical :: down_flux
1921 if (continuous)
then
1922 call interface_scalar(nk, hl, tl, til, 2, h_neglect)
1923 call interface_scalar(nk, hr, tr, tir, 2, h_neglect)
1924 call ppm_left_right_edge_values(nk, tl, til, al_l, ar_l)
1925 call ppm_left_right_edge_values(nk, tr, tir, al_r, ar_r)
1927 ppoly_r_coeffs_l(:,:) = 0.
1928 ppoly_r_coeffs_r(:,:) = 0.
1932 call build_reconstructions_1d( remap_cs, nk, hl, tl, ppoly_r_coeffs_l, tid_l, &
1933 ppoly_r_s_l, imethod, h_neglect, h_neglect_edge )
1934 call build_reconstructions_1d( remap_cs, nk, hr, tr, ppoly_r_coeffs_r, tid_r, &
1935 ppoly_r_s_r, imethod, h_neglect, h_neglect_edge )
1938 do k_sublayer = 1, nsurf-1
1939 if (heff(k_sublayer) == 0.)
then
1940 flx(k_sublayer) = 0.
1942 if (continuous)
then
1943 klb = kol(k_sublayer+1)
1944 t_left_bottom = ( 1. - pil(k_sublayer+1) ) * til(klb) + pil(k_sublayer+1) * til(klb+1)
1945 klt = kol(k_sublayer)
1946 t_left_top = ( 1. - pil(k_sublayer) ) * til(klt) + pil(k_sublayer) * til(klt+1)
1947 t_left_layer = ppm_ave(pil(k_sublayer), pil(k_sublayer+1) + real(klb-klt), &
1948 al_l(klt), ar_l(klt), tl(klt))
1950 krb = kor(k_sublayer+1)
1951 t_right_bottom = ( 1. - pir(k_sublayer+1) ) * tir(krb) + pir(k_sublayer+1) * tir(krb+1)
1952 krt = kor(k_sublayer)
1953 t_right_top = ( 1. - pir(k_sublayer) ) * tir(krt) + pir(k_sublayer) * tir(krt+1)
1954 t_right_layer = ppm_ave(pir(k_sublayer), pir(k_sublayer+1) + real(krb-krt), &
1955 al_r(krt), ar_r(krt), tr(krt))
1956 dt_top = t_right_top - t_left_top
1957 dt_bottom = t_right_bottom - t_left_bottom
1958 dt_ave = 0.5 * ( dt_top + dt_bottom )
1959 dt_layer = t_right_layer - t_left_layer
1960 if (signum(1.,dt_top) * signum(1.,dt_bottom) <= 0. .or. signum(1.,dt_ave) * signum(1.,dt_layer) <= 0.)
then
1965 flx(k_sublayer) = dt_ave * heff(k_sublayer)
1968 call neutral_surface_t_eval(nk, nsurf, k_sublayer, kol, pil, tl, tid_l, deg, imethod, &
1969 ppoly_r_coeffs_l, t_left_top, t_left_bottom, t_left_sub, &
1970 t_left_top_int, t_left_bot_int, t_left_layer)
1971 call neutral_surface_t_eval(nk, nsurf, k_sublayer, kor, pir, tr, tid_r, deg, imethod, &
1972 ppoly_r_coeffs_r, t_right_top, t_right_bottom, t_right_sub, &
1973 t_right_top_int, t_right_bot_int, t_right_layer)
1975 dt_top = t_right_top - t_left_top
1976 dt_bottom = t_right_bottom - t_left_bottom
1977 dt_sublayer = t_right_sub - t_left_sub
1978 dt_top_int = t_right_top_int - t_left_top_int
1979 dt_bot_int = t_right_bot_int - t_left_bot_int
1983 down_flux = dt_top <= 0. .and. dt_bottom <= 0. .and. &
1984 dt_sublayer <= 0. .and. dt_top_int <= 0. .and. &
1986 down_flux = down_flux .or. &
1987 (dt_top >= 0. .and. dt_bottom >= 0. .and. &
1988 dt_sublayer >= 0. .and. dt_top_int >= 0. .and. &
1991 flx(k_sublayer) = dt_sublayer * heff(k_sublayer)
1993 flx(k_sublayer) = 0.
1999 end subroutine neutral_surface_flux
2002 subroutine neutral_surface_t_eval(nk, ns, k_sub, Ks, Ps, T_mean, T_int, deg, iMethod, T_poly, &
2003 T_top, T_bot, T_sub, T_top_int, T_bot_int, T_layer)
2004 integer,
intent(in ) :: nk
2005 integer,
intent(in ) :: ns
2006 integer,
intent(in ) :: k_sub
2007 integer,
dimension(ns),
intent(in ) :: Ks
2008 real,
dimension(ns),
intent(in ) :: Ps
2009 real,
dimension(nk),
intent(in ) :: T_mean
2010 real,
dimension(nk,2),
intent(in ) :: T_int
2011 integer,
intent(in ) :: deg
2012 integer,
intent(in ) :: iMethod
2013 real,
dimension(nk,deg+1),
intent(in ) :: T_poly
2014 real,
intent( out) :: T_top
2015 real,
intent( out) :: T_bot
2016 real,
intent( out) :: T_sub
2017 real,
intent( out) :: T_top_int
2018 real,
intent( out) :: T_bot_int
2019 real,
intent( out) :: T_layer
2021 integer :: kl, ks_top, ks_bot
2025 if ( ks(ks_top) /= ks(ks_bot) )
then
2026 call mom_error(fatal,
"Neutral surfaces span more than one layer")
2030 if ( (ps(ks_top) == 0.) .and. (ps(ks_bot) == 1.))
then
2035 if ( (kl > 1) .and. (ps(ks_top) == 0.) )
then
2036 t_top = t_int(kl-1,2)
2038 t_top = evaluation_polynomial( t_poly(kl,:), deg+1, ps(ks_top) )
2041 if ( (kl < nk) .and. (ps(ks_bot) == 1.) )
then
2042 t_bot = t_int(kl+1,1)
2044 t_bot = evaluation_polynomial( t_poly(kl,:), deg+1, ps(ks_bot) )
2047 t_sub = average_value_ppoly(nk, t_mean, t_int, t_poly, imethod, kl, ps(ks_top), ps(ks_bot))
2048 t_top_int = evaluation_polynomial( t_poly(kl,:), deg+1, ps(ks_top))
2049 t_bot_int = evaluation_polynomial( t_poly(kl,:), deg+1, ps(ks_bot))
2050 t_layer = t_mean(kl)
2052 end subroutine neutral_surface_t_eval
2055 subroutine ppm_left_right_edge_values(nk, Tl, Ti, aL, aR)
2056 integer,
intent(in) :: nk
2057 real,
dimension(nk),
intent(in) :: Tl
2058 real,
dimension(nk+1),
intent(in) :: Ti
2059 real,
dimension(nk),
intent(inout) :: aL
2060 real,
dimension(nk),
intent(inout) :: aR
2067 if ( signum(1., ar(k) - tl(k))*signum(1., tl(k) - al(k)) <= 0.0 )
then
2070 elseif ( sign(3., ar(k) - al(k)) * ( (tl(k) - al(k)) + (tl(k) - ar(k))) > abs(ar(k) - al(k)) )
then
2071 al(k) = tl(k) + 2.0 * ( tl(k) - ar(k) )
2072 elseif ( sign(3., ar(k) - al(k)) * ( (tl(k) - al(k)) + (tl(k) - ar(k))) < -abs(ar(k) - al(k)) )
then
2073 ar(k) = tl(k) + 2.0 * ( tl(k) - al(k) )
2076 end subroutine ppm_left_right_edge_values
2079 logical function neutral_diffusion_unit_tests(verbose)
2080 logical,
intent(in) :: verbose
2082 neutral_diffusion_unit_tests = .false. .or. &
2083 ndiff_unit_tests_continuous(verbose) .or. ndiff_unit_tests_discontinuous(verbose)
2085 end function neutral_diffusion_unit_tests
2088 logical function ndiff_unit_tests_continuous(verbose)
2089 logical,
intent(in) :: verbose
2091 integer,
parameter :: nk = 4
2092 real,
dimension(nk+1) :: til, tir1, tir2, tir4, tio
2093 real,
dimension(nk) :: tl
2094 real,
dimension(nk+1) :: sil
2095 real,
dimension(nk+1) :: pil, pir4
2096 real,
dimension(2*nk+2) :: pilro, pirlo
2097 integer,
dimension(2*nk+2) :: kol, kor
2098 real,
dimension(2*nk+1) :: heff
2099 real,
dimension(2*nk+1) :: flx
2108 ndiff_unit_tests_continuous = .false.
2109 write(stdout,*)
'==== MOM_neutral_diffusion: ndiff_unit_tests_continuous ='
2111 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2112 test_fv_diff(v,1.,1.,1., 0.,1.,2., 1.,
'FV: Straight line on uniform grid')
2113 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2114 test_fv_diff(v,1.,1.,0., 0.,4.,8., 7.,
'FV: Vanished right cell')
2115 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2116 test_fv_diff(v,0.,1.,1., 0.,4.,8., 7.,
'FV: Vanished left cell')
2117 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2118 test_fv_diff(v,1.,2.,4., 0.,3.,9., 4.,
'FV: Stretched grid')
2119 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2120 test_fv_diff(v,2.,0.,2., 0.,1.,2., 0.,
'FV: Vanished middle cell')
2121 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2122 test_fv_diff(v,0.,1.,0., 0.,1.,2., 2.,
'FV: Vanished on both sides')
2123 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2124 test_fv_diff(v,1.,0.,0., 0.,1.,2., 0.,
'FV: Two vanished cell sides')
2125 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2126 test_fv_diff(v,0.,0.,0., 0.,1.,2., 0.,
'FV: All vanished cells')
2128 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2129 test_fvlsq_slope(v,1.,1.,1., 0.,1.,2., 1.,
'LSQ: Straight line on uniform grid')
2130 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2131 test_fvlsq_slope(v,1.,1.,0., 0.,1.,2., 1.,
'LSQ: Vanished right cell')
2132 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2133 test_fvlsq_slope(v,0.,1.,1., 0.,1.,2., 1.,
'LSQ: Vanished left cell')
2134 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2135 test_fvlsq_slope(v,1.,2.,4., 0.,3.,9., 2.,
'LSQ: Stretched grid')
2136 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2137 test_fvlsq_slope(v,1.,0.,1., 0.,1.,2., 2.,
'LSQ: Vanished middle cell')
2138 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2139 test_fvlsq_slope(v,0.,1.,0., 0.,1.,2., 0.,
'LSQ: Vanished on both sides')
2140 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2141 test_fvlsq_slope(v,1.,0.,0., 0.,1.,2., 0.,
'LSQ: Two vanished cell sides')
2142 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2143 test_fvlsq_slope(v,0.,0.,0., 0.,1.,2., 0.,
'LSQ: All vanished cells')
2145 call interface_scalar(4, (/10.,10.,10.,10./), (/24.,18.,12.,6./), tio, 1, h_neglect)
2148 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2149 test_data1d(v,5, tio, (/24.,22.5,15.,7.5,6./),
'Linear profile, linear interface temperatures')
2150 call interface_scalar(4, (/10.,10.,10.,10./), (/24.,18.,12.,6./), tio, 2, h_neglect)
2151 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2152 test_data1d(v,5, tio, (/24.,22.,15.,8.,6./),
'Linear profile, PPM interface temperatures')
2154 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2155 test_ifndp(v,-1.0, 0., 1.0, 1.0, 0.5,
'Check mid-point')
2156 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2157 test_ifndp(v, 0.0, 0., 1.0, 1.0, 0.0,
'Check bottom')
2158 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2159 test_ifndp(v, 0.1, 0., 1.1, 1.0, 0.0,
'Check below')
2160 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2161 test_ifndp(v,-1.0, 0., 0.0, 1.0, 1.0,
'Check top')
2162 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2163 test_ifndp(v,-1.0, 0., -0.1, 1.0, 1.0,
'Check above')
2164 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2165 test_ifndp(v,-1.0, 0., 3.0, 1.0, 0.25,
'Check 1/4')
2166 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2167 test_ifndp(v,-3.0, 0., 1.0, 1.0, 0.75,
'Check 3/4')
2168 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2169 test_ifndp(v, 1.0, 0., 1.0, 1.0, 0.0,
'Check dRho=0 below')
2170 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2171 test_ifndp(v,-1.0, 0., -1.0, 1.0, 1.0,
'Check dRho=0 above')
2172 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2173 test_ifndp(v, 0.0, 0., 0.0, 1.0, 0.5,
'Check dRho=0 mid')
2174 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. &
2175 test_ifndp(v,-2.0, .5, 5.0, 0.5, 0.5,
'Check dP=0')
2178 call find_neutral_surface_positions_continuous(3, &
2179 (/0.,10.,20.,30./), (/22.,18.,14.,10./), (/0.,0.,0.,0./), &
2180 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2181 (/0.,10.,20.,30./), (/22.,18.,14.,10./), (/0.,0.,0.,0./), &
2182 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2183 pilro, pirlo, kol, kor, heff)
2184 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_nsp(v, 8, kol, kor, pilro, pirlo, heff, &
2185 (/1,1,2,2,3,3,3,3/), &
2186 (/1,1,2,2,3,3,3,3/), &
2187 (/0.,0.,0.,0.,0.,0.,1.,1./), &
2188 (/0.,0.,0.,0.,0.,0.,1.,1./), &
2189 (/0.,10.,0.,10.,0.,10.,0./), &
2190 'Identical columns')
2191 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_data1d(v, 8, &
2192 absolute_positions(3, 8, (/0.,10.,20.,30./), kol, pilro), &
2193 (/0.,0.,10.,10.,20.,20.,30.,30./),
'... left positions')
2194 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_data1d(v, 8, &
2195 absolute_positions(3, 8, (/0.,10.,20.,30./), kor, pirlo), &
2196 (/0.,0.,10.,10.,20.,20.,30.,30./),
'... right positions')
2197 call neutral_surface_flux(3, 2*3+2, 2, (/10.,10.,10./), (/10.,10.,10./), &
2198 (/20.,16.,12./), (/20.,16.,12./), &
2199 pilro, pirlo, kol, kor, heff, flx, .true., h_neglect)
2200 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_data1d(v, 7, flx, &
2201 (/0.,0.,0.,0.,0.,0.,0./),
'Identical columns, rho flux (=0)')
2202 call neutral_surface_flux(3, 2*3+2, 2, (/10.,10.,10./), (/10.,10.,10./), &
2203 (/-1.,-1.,-1./), (/1.,1.,1./), &
2204 pilro, pirlo, kol, kor, heff, flx, .true., h_neglect)
2205 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_data1d(v, 7, flx, &
2206 (/0.,20.,0.,20.,0.,20.,0./),
'Identical columns, S flux')
2209 call find_neutral_surface_positions_continuous(3, &
2210 (/0.,10.,20.,30./), (/22.,18.,14.,10./), (/0.,0.,0.,0./), &
2211 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2212 (/0.,10.,20.,30./), (/20.,16.,12.,8./), (/0.,0.,0.,0./), &
2213 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2214 pilro, pirlo, kol, kor, heff)
2215 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_nsp(v, 8, kol, kor, pilro, pirlo, heff, &
2216 (/1,1,2,2,3,3,3,3/), &
2217 (/1,1,1,2,2,3,3,3/), &
2218 (/0.,0.5,0.,0.5,0.,0.5,1.,1./), &
2219 (/0.,0.,0.5,0.,0.5,0.,0.5,1./), &
2220 (/0.,5.,5.,5.,5.,5.,0./), &
2221 'Right column slightly cooler')
2222 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_data1d(v, 8, &
2223 absolute_positions(3, 8, (/0.,10.,20.,30./), kol, pilro), &
2224 (/0.,5.,10.,15.,20.,25.,30.,30./),
'... left positions')
2225 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_data1d(v, 8, &
2226 absolute_positions(3, 8, (/0.,10.,20.,30./), kor, pirlo), &
2227 (/0.,0.,5.,10.,15.,20.,25.,30./),
'... right positions')
2230 call find_neutral_surface_positions_continuous(3, &
2231 (/0.,10.,20.,30./), (/22.,18.,14.,10./), (/0.,0.,0.,0./), &
2232 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2233 (/0.,10.,20.,30./), (/24.,20.,16.,12./), (/0.,0.,0.,0./), &
2234 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2235 pilro, pirlo, kol, kor, heff)
2236 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_nsp(v, 8, kol, kor, pilro, pirlo, heff, &
2237 (/1,1,1,2,2,3,3,3/), &
2238 (/1,1,2,2,3,3,3,3/), &
2239 (/0.,0.,0.5,0.,0.5,0.,0.5,1./), &
2240 (/0.,0.5,0.,0.5,0.,0.5,1.,1./), &
2241 (/0.,5.,5.,5.,5.,5.,0./), &
2242 'Right column slightly warmer')
2245 call find_neutral_surface_positions_continuous(3, &
2246 (/0.,10.,20.,30./), (/22.,18.,14.,10./), (/0.,0.,0.,0./), &
2247 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2248 (/0.,10.,20.,30./), (/16.,12.,8.,4./), (/0.,0.,0.,0./), &
2249 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2250 pilro, pirlo, kol, kor, heff)
2251 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_nsp(v, 8, kol, kor, pilro, pirlo, heff, &
2252 (/1,2,2,3,3,3,3,3/), &
2253 (/1,1,1,1,2,2,3,3/), &
2254 (/0.,0.,0.5,0.,0.5,1.,1.,1./), &
2255 (/0.,0.,0.,0.5,0.,0.5,0.,1./), &
2256 (/0.,0.,5.,5.,5.,0.,0./), &
2257 'Right column somewhat cooler')
2260 call find_neutral_surface_positions_continuous(3, &
2261 (/0.,10.,20.,30./), (/22.,18.,14.,10./), (/0.,0.,0.,0./), &
2262 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2263 (/0.,10.,20.,30./), (/9.,7.,5.,3./), (/0.,0.,0.,0./), &
2264 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2265 pilro, pirlo, kol, kor, heff)
2266 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_nsp(v, 8, kol, kor, pilro, pirlo, heff, &
2267 (/1,2,3,3,3,3,3,3/), &
2268 (/1,1,1,1,1,2,3,3/), &
2269 (/0.,0.,0.,1.,1.,1.,1.,1./), &
2270 (/0.,0.,0.,0.,0.,0.,0.,1./), &
2271 (/0.,0.,0.,0.,0.,0.,0./), &
2272 'Right column much cooler')
2275 call find_neutral_surface_positions_continuous(3, &
2276 (/0.,10.,20.,30./), (/22.,18.,14.,10./), (/0.,0.,0.,0./), &
2277 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2278 (/0.,10.,20.,30./), (/14.,14.,10.,2./), (/0.,0.,0.,0./), &
2279 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2280 pilro, pirlo, kol, kor, heff)
2281 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_nsp(v, 8, kol, kor, pilro, pirlo, heff, &
2282 (/1,2,3,3,3,3,3,3/), &
2283 (/1,1,1,1,2,3,3,3/), &
2284 (/0.,0.,0.,0.,0.,1.,1.,1./), &
2285 (/0.,0.,0.,0.,0.,0.,0.,1./), &
2286 (/0.,0.,0.,0.,10.,0.,0./), &
2287 'Right column with mixed layer')
2290 call find_neutral_surface_positions_continuous(3, &
2291 (/0.,10.,20.,30./), (/14.,14.,10.,2./), (/0.,0.,0.,0./), &
2292 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2293 (/0.,10.,20.,30./), (/14.,14.,10.,2./), (/0.,0.,0.,0./), &
2294 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2295 pilro, pirlo, kol, kor, heff)
2296 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_nsp(v, 8, kol, kor, pilro, pirlo, heff, &
2297 (/1,1,2,2,3,3,3,3/), &
2298 (/1,1,2,2,3,3,3,3/), &
2299 (/0.,0.,0.,0.,0.,0.,1.,1./), &
2300 (/0.,0.,0.,0.,0.,0.,1.,1./), &
2301 (/0.,10.,0.,10.,0.,10.,0./), &
2302 'Identical columns with mixed layer')
2305 call find_neutral_surface_positions_continuous(3, &
2306 (/0.,10.,20.,30./), (/14.,14.,10.,2./), (/0.,0.,0.,0./), &
2307 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2308 (/0.,10.,20.,30./), (/10.,14.,12.,4./), (/0.,0.,0.,0./), &
2309 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2310 pilro, pirlo, kol, kor, heff)
2311 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_nsp(v, 8, kol, kor, pilro, pirlo, heff, &
2312 (/1,2,3,3,3,3,3,3/), &
2313 (/1,1,1,2,3,3,3,3/), &
2314 (/0.,0.,0.,0.,0.,0.,.75,1./), &
2315 (/0.,0.,0.,0.,0.,0.25,1.,1./), &
2316 (/0.,0.,0.,0.,0.,7.5,0./), &
2317 'Right column with unstable mixed layer')
2320 call find_neutral_surface_positions_continuous(3, &
2321 (/0.,10.,20.,30./), (/10.,14.,12.,4./), (/0.,0.,0.,0./), &
2322 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2323 (/0.,10.,20.,30./), (/14.,14.,10.,2./), (/0.,0.,0.,0./), &
2324 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2325 pilro, pirlo, kol, kor, heff)
2326 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_nsp(v, 8, kol, kor, pilro, pirlo, heff, &
2327 (/1,1,1,2,3,3,3,3/), &
2328 (/1,2,3,3,3,3,3,3/), &
2329 (/0.,0.,0.,0.,0.,0.25,1.,1./), &
2330 (/0.,0.,0.,0.,0.,0.,.75,1./), &
2331 (/0.,0.,0.,0.,0.,7.5,0./), &
2332 'Left column with unstable mixed layer')
2335 call find_neutral_surface_positions_continuous(3, &
2336 (/0.,10.,20.,30./), (/8.,12.,10.,2./), (/0.,0.,0.,0./), &
2337 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2338 (/0.,10.,20.,30./), (/10.,14.,12.,4./), (/0.,0.,0.,0./), &
2339 (/-1.,-1.,-1.,-1./), (/1.,1.,1.,1./), &
2340 pilro, pirlo, kol, kor, heff)
2341 ndiff_unit_tests_continuous = ndiff_unit_tests_continuous .or. test_nsp(v, 8, kol, kor, pilro, pirlo, heff, &
2342 (/1,1,1,1,2,3,3,3/), &
2343 (/1,2,3,3,3,3,3,3/), &
2344 (/0.,0.,0.,0.,0.,0.,0.75,1./), &
2345 (/0.,0.,0.,0.5,0.5,0.5,1.,1./), &
2346 (/0.,0.,0.,0.,0.,6.,0./), &
2347 'Two unstable mixed layers')
2349 if (.not. ndiff_unit_tests_continuous)
write(stdout,*)
'Pass'
2351 end function ndiff_unit_tests_continuous
2353 logical function ndiff_unit_tests_discontinuous(verbose)
2354 logical,
intent(in) :: verbose
2356 integer,
parameter :: nk = 3
2357 integer,
parameter :: ns = nk*4
2358 real,
dimension(nk) :: sl, sr, tl, tr
2359 real,
dimension(nk) :: hl, hr
2360 real,
dimension(nk,2) :: til, sil, tir, sir
2361 real,
dimension(nk,2) :: pres_l, pres_r
2362 integer,
dimension(ns) :: kol, kor
2363 real,
dimension(ns) :: pol, por
2364 real,
dimension(ns-1) :: heff, flx
2368 real,
dimension(nk,2) :: ppoly_t_l, ppoly_t_r
2369 real,
dimension(nk,2) :: ppoly_s_l, ppoly_s_r
2370 real,
dimension(nk,2) :: drdt
2372 real,
dimension(nk,2) :: drds
2374 logical,
dimension(nk) :: stable_l, stable_r
2376 integer :: ns_l, ns_r
2381 ndiff_unit_tests_discontinuous = .false.
2382 write(stdout,*)
'==== MOM_neutral_diffusion: ndiff_unit_tests_discontinuous ='
2387 call eos_manual_init(cs%EOS, form_of_eos=eos_linear, drho_dt=-1., drho_ds=0.)
2388 sl(:) = 0. ; sr(:) = 0. ; ; sil(:,:) = 0. ; sir(:,:) = 0.
2389 ppoly_t_l(:,:) = 0.; ppoly_t_r(:,:) = 0.
2390 ppoly_s_l(:,:) = 0.; ppoly_s_r(:,:) = 0.
2394 hl = (/10.,10.,10./) ; hr = (/10.,10.,10./)
2395 pres_l(1,1) = 0. ; pres_l(1,2) = hl(1) ; pres_r(1,1) = 0. ; pres_r(1,2) = hr(1)
2397 pres_l(k,1) = pres_l(k-1,2)
2398 pres_l(k,2) = pres_l(k,1) + hl(k)
2399 pres_r(k,1) = pres_r(k-1,2)
2400 pres_r(k,2) = pres_r(k,1) + hr(k)
2402 cs%delta_rho_form =
'mid_pressure'
2403 cs%neutral_pos_method = 1
2405 til(1,:) = (/ 22.00, 18.00 /); til(2,:) = (/ 18.00, 14.00 /); til(3,:) = (/ 14.00, 10.00 /);
2406 tir(1,:) = (/ 22.00, 18.00 /); tir(2,:) = (/ 18.00, 14.00 /); tir(3,:) = (/ 14.00, 10.00 /);
2407 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2408 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2409 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2410 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2411 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2412 (/ 1, 1, 1, 1, 2, 2, 2, 2, 3, 3, 3, 3 /), &
2413 (/ 1, 1, 1, 2, 2, 2, 2, 3, 3, 3, 3, 3 /), &
2414 (/ 0.00, 0.00, 1.00, 1.00, 0.00, 0.00, 1.00, 1.00, 0.00, 0.00, 1.00, 1.00 /), &
2415 (/ 0.00, 0.00, 1.00, 0.00, 0.00, 0.00, 1.00, 0.00, 0.00, 0.00, 1.00, 1.00 /), &
2416 (/ 0.00, 10.00, 0.00, 0.00, 0.00, 10.00, 0.00, 0.00, 0.00, 10.00, 0.00 /), &
2417 'Identical Columns')
2419 til(1,:) = (/ 22.00, 18.00 /); til(2,:) = (/ 18.00, 14.00 /); til(3,:) = (/ 14.00, 10.00 /);
2420 tir(1,:) = (/ 20.00, 16.00 /); tir(2,:) = (/ 16.00, 12.00 /); tir(3,:) = (/ 12.00, 8.00 /);
2421 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2422 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2423 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2424 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2425 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2426 (/ 1, 1, 1, 2, 2, 2, 2, 3, 3, 3, 3, 3 /), &
2427 (/ 1, 1, 1, 1, 1, 2, 2, 2, 2, 3, 3, 3 /), &
2428 (/ 0.00, 0.50, 1.00, 0.00, 0.50, 0.50, 1.00, 0.00, 0.50, 0.50, 1.00, 1.00 /), &
2429 (/ 0.00, 0.00, 0.50, 0.50, 1.00, 0.00, 0.50, 0.50, 1.00, 0.00, 0.50, 1.00 /), &
2430 (/ 0.00, 5.00, 0.00, 5.00, 0.00, 5.00, 0.00, 5.00, 0.00, 5.00, 0.00 /), &
2431 'Right slightly cooler')
2433 til(1,:) = (/ 20.00, 16.00 /); til(2,:) = (/ 16.00, 12.00 /); til(3,:) = (/ 12.00, 8.00 /);
2434 tir(1,:) = (/ 22.00, 18.00 /); tir(2,:) = (/ 18.00, 14.00 /); tir(3,:) = (/ 14.00, 10.00 /);
2435 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2436 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2437 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2438 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2439 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2440 (/ 1, 1, 1, 1, 1, 2, 2, 2, 2, 3, 3, 3 /), &
2441 (/ 1, 1, 1, 2, 2, 2, 2, 3, 3, 3, 3, 3 /), &
2442 (/ 0.00, 0.00, 0.50, 0.50, 1.00, 0.00, 0.50, 0.50, 1.00, 0.00, 0.50, 1.00 /), &
2443 (/ 0.00, 0.50, 1.00, 0.00, 0.50, 0.50, 1.00, 0.00, 0.50, 0.50, 1.00, 1.00 /), &
2444 (/ 0.00, 5.00, 0.00, 5.00, 0.00, 5.00, 0.00, 5.00, 0.00, 5.00, 0.00 /), &
2445 'Left slightly cooler')
2447 til(1,:) = (/ 22.00, 20.00 /); til(2,:) = (/ 18.00, 16.00 /); til(3,:) = (/ 14.00, 12.00 /);
2448 tir(1,:) = (/ 32.00, 24.00 /); tir(2,:) = (/ 22.00, 14.00 /); tir(3,:) = (/ 12.00, 4.00 /);
2449 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2450 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2451 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2452 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2453 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2454 (/ 1, 1, 1, 1, 1, 2, 2, 3, 3, 3, 3, 3 /), &
2455 (/ 1, 1, 2, 2, 2, 2, 2, 2, 3, 3, 3, 3 /), &
2456 (/ 0.00, 0.00, 0.00, 0.00, 1.00, 0.00, 1.00, 0.00, 0.00, 1.00, 1.00, 1.00 /), &
2457 (/ 0.00, 1.00, 0.00, 0.00, 0.25, 0.50, 0.75, 1.00, 0.00, 0.00, 0.00, 1.00 /), &
2458 (/ 0.00, 0.00, 0.00, 4.00, 0.00, 4.00, 0.00, 0.00, 0.00, 0.00, 0.00 /), &
2459 'Right more strongly stratified')
2461 til(1,:) = (/ 22.00, 18.00 /); til(2,:) = (/ 18.00, 14.00 /); til(3,:) = (/ 14.00, 10.00 /);
2462 tir(1,:) = (/ 14.00, 14.00 /); tir(2,:) = (/ 14.00, 14.00 /); tir(3,:) = (/ 12.00, 8.00 /);
2463 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2464 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2465 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2466 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2467 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2468 (/ 1, 1, 1, 1, 1, 1, 2, 2, 3, 3, 3, 3 /), &
2469 (/ 1, 1, 2, 2, 3, 3, 3, 3, 3, 3, 3, 3 /), &
2470 (/ 0.00, 0.00, 0.00, 0.00, 0.00, 1.00, 0.00, 1.00, 0.00, 0.50, 1.00, 1.00 /), &
2471 (/ 0.00, 1.00, 0.00, 1.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.50, 1.00 /), &
2472 (/ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 5.00, 0.00 /), &
2473 'Deep Mixed layer on the right')
2475 til(1,:) = (/ 14.00, 14.00 /); til(2,:) = (/ 14.00, 12.00 /); til(3,:) = (/ 10.00, 8.00 /);
2476 tir(1,:) = (/ 14.00, 14.00 /); tir(2,:) = (/ 14.00, 14.00 /); tir(3,:) = (/ 14.00, 14.00 /);
2477 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2478 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2479 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2480 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2481 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2482 (/ 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 3 /), &
2483 (/ 1, 1, 1, 1, 2, 2, 3, 3, 3, 3, 3, 3 /), &
2484 (/ 0.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00 /), &
2485 (/ 0.00, 0.00, 0.00, 1.00, 0.00, 1.00, 0.00, 1.00, 1.00, 1.00, 1.00, 1.00 /), &
2486 (/ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00 /), &
2487 'Right unstratified column')
2489 til(1,:) = (/ 14.00, 14.00 /); til(2,:) = (/ 14.00, 12.00 /); til(3,:) = (/ 10.00, 8.00 /);
2490 tir(1,:) = (/ 14.00, 14.00 /); tir(2,:) = (/ 14.00, 14.00 /); tir(3,:) = (/ 12.00, 4.00 /);
2491 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2492 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2493 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2494 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2495 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2496 (/ 1, 1, 1, 1, 1, 1, 2, 2, 2, 3, 3, 3 /), &
2497 (/ 1, 1, 1, 1, 2, 2, 3, 3, 3, 3, 3, 3 /), &
2498 (/ 0.00, 1.00, 1.00, 1.00, 1.00, 1.00, 0.00, 1.00, 1.00, 0.00, 1.00, 1.00 /), &
2499 (/ 0.00, 0.00, 0.00, 1.00, 0.00, 1.00, 0.00, 0.00, 0.00, 0.25, 0.50, 1.00 /), &
2500 (/ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 4.00, 0.00 /), &
2501 'Right unstratified column')
2503 til(1,:) = (/ 14.00, 14.00 /); til(2,:) = (/ 14.00, 10.00 /); til(3,:) = (/ 10.00, 2.00 /);
2504 tir(1,:) = (/ 14.00, 14.00 /); tir(2,:) = (/ 14.00, 10.00 /); tir(3,:) = (/ 10.00, 2.00 /);
2505 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2506 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2507 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2508 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2509 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2510 (/ 1, 1, 1, 1, 2, 2, 2, 2, 3, 3, 3, 3 /), &
2511 (/ 1, 1, 1, 1, 2, 2, 2, 3, 3, 3, 3, 3 /), &
2512 (/ 0.00, 1.00, 1.00, 1.00, 0.00, 0.00, 1.00, 1.00, 0.00, 0.00, 1.00, 1.00 /), &
2513 (/ 0.00, 0.00, 0.00, 1.00, 0.00, 0.00, 1.00, 0.00, 0.00, 0.00, 1.00, 1.00 /), &
2514 (/ 0.00, 0.00, 0.00, 0.00, 0.00, 10.00, 0.00, 0.00, 0.00, 10.00, 0.00 /), &
2515 'Identical columns with mixed layer')
2517 til(1,:) = (/ 14.00, 12.00 /); til(2,:) = (/ 10.00, 10.00 /); til(3,:) = (/ 8.00, 2.00 /);
2518 tir(1,:) = (/ 14.00, 12.00 /); tir(2,:) = (/ 12.00, 8.00 /); tir(3,:) = (/ 8.00, 2.00 /);
2519 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2520 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2521 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2522 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2523 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2524 (/ 1, 1, 1, 1, 2, 2, 3, 3, 3, 3, 3, 3 /), &
2525 (/ 1, 1, 1, 2, 2, 2, 2, 2, 2, 3, 3, 3 /), &
2526 (/ 0.00, 0.00, 1.00, 1.00, 0.00, 1.00, 0.00, 0.00, 0.00, 0.00, 1.00, 1.00 /), &
2527 (/ 0.00, 0.00, 1.00, 0.00, 0.00, 0.00, 0.00, 1.00, 1.00, 0.00, 1.00, 1.00 /), &
2528 (/ 0.00, 10.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 10.00, 0.00 /), &
2529 'Left interior unstratified')
2531 til(1,:) = (/ 12.00, 12.00 /); til(2,:) = (/ 12.00, 10.00 /); til(3,:) = (/ 10.00, 6.00 /);
2532 tir(1,:) = (/ 12.00, 10.00 /); tir(2,:) = (/ 10.00, 12.00 /); tir(3,:) = (/ 8.00, 4.00 /);
2533 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2534 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2535 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2536 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2537 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2538 (/ 1, 1, 2, 2, 2, 3, 3, 3, 3, 3, 3, 3 /), &
2539 (/ 1, 1, 1, 1, 1, 1, 2, 2, 3, 3, 3, 3 /), &
2540 (/ 0.00, 1.00, 0.00, 0.00, 1.00, 0.00, 0.00, 0.00, 0.00, 0.50, 1.00, 1.00 /), &
2541 (/ 0.00, 0.00, 0.00, 0.00, 1.00, 1.00, 0.00, 1.00, 0.00, 0.00, 0.50, 1.00 /), &
2542 (/ 0.00, 0.00, 0.00, 10.00, 0.00, 0.00, 0.00, 0.00, 0.00, 5.00, 0.00 /), &
2543 'Left mixed layer, Right unstable interior')
2545 til(1,:) = (/ 14.00, 14.00 /); til(2,:) = (/ 10.00, 10.00 /); til(3,:) = (/ 8.00, 6.00 /);
2546 tir(1,:) = (/ 10.00, 14.00 /); tir(2,:) = (/ 16.00, 16.00 /); tir(3,:) = (/ 12.00, 4.00 /);
2547 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2548 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2549 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2550 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2551 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2552 (/ 1, 1, 2, 2, 2, 2, 2, 2, 3, 3, 3, 3 /), &
2553 (/ 1, 1, 1, 1, 1, 1, 2, 2, 3, 3, 3, 3 /), &
2554 (/ 0.00, 1.00, 0.00, 1.00, 1.00, 1.00, 1.00, 1.00, 0.00, 0.00, 1.00, 1.00 /), &
2555 (/ 0.00, 0.00, 0.00, 0.00, 0.00, 1.00, 0.00, 1.00, 0.00, 0.50, 0.75, 1.00 /), &
2556 (/ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 4.00, 0.00 /), &
2557 'Left thick mixed layer, Right unstable mixed')
2559 til(1,:) = (/ 8.00, 12.00 /); til(2,:) = (/ 12.00, 10.00 /); til(3,:) = (/ 8.00, 4.00 /);
2560 tir(1,:) = (/ 10.00, 14.00 /); tir(2,:) = (/ 14.00, 12.00 /); tir(3,:) = (/ 10.00, 6.00 /);
2561 call mark_unstable_cells( cs, nk, til, sil, pres_l, stable_l )
2562 call mark_unstable_cells( cs, nk, tir, sir, pres_r, stable_r )
2563 call find_neutral_surface_positions_discontinuous(cs, nk, pres_l, hl, til, sil, ppoly_t_l, ppoly_s_l, stable_l, &
2564 pres_r, hr, tir, sir, ppoly_t_r, ppoly_s_r, stable_r, pol, por, kol, kor, heff)
2565 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. test_nsp(v, 12, kol, kor, pol, por, heff, &
2566 (/ 1, 1, 1, 1, 2, 2, 2, 2, 3, 3, 3, 3 /), &
2567 (/ 1, 1, 1, 1, 2, 2, 2, 3, 3, 3, 3, 3 /), &
2568 (/ 0.00, 1.00, 1.00, 1.00, 0.00, 0.00, 0.00, 1.00, 0.00, 0.00, 0.50, 1.00 /), &
2569 (/ 0.00, 0.00, 0.00, 1.00, 0.00, 1.00, 1.00, 0.00, 0.00, 0.50, 1.00, 1.00 /), &
2570 (/ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 5.00, 0.00 /), &
2571 'Unstable mixed layers, left cooler')
2573 call eos_manual_init(cs%EOS, form_of_eos = eos_linear, drho_dt = -1., drho_ds = 2.)
2580 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. (test_rnp(0.5, &
2581 find_neutral_pos_linear(cs, 0., 10., 35., -0.2, 0., &
2582 -0.2, 0., -0.2, 0., &
2583 (/12.,-4./), (/34.,0./)),
"Temp Uniform Linearized Alpha/Beta"))
2585 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. (test_rnp(0.5, &
2586 find_neutral_pos_linear(cs, 0., 10., 35., 0., 0.8, &
2588 (/12.,0./), (/34.,2./)),
"Salt Uniform Linearized Alpha/Beta"))
2590 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. (test_rnp(0.5, &
2591 find_neutral_pos_linear(cs, 0., 10., 35., -0.5, 0.5, &
2592 -0.5, 0.5, -0.5, 0.5, &
2593 (/12.,-4./), (/34.,2./)),
"Temp/salt Uniform Linearized Alpha/Beta"))
2595 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. (test_rnp(0.5, &
2596 find_neutral_pos_linear(cs, 0., 10., 35., -0.2, 0., &
2597 -0.4, 0., -0.6, 0., &
2598 (/12.,-4./), (/34.,0./)),
"Temp stratified Linearized Alpha/Beta"))
2600 ndiff_unit_tests_discontinuous = ndiff_unit_tests_discontinuous .or. (test_rnp(0.5, &
2601 find_neutral_pos_linear(cs, 0., 10., 35., 0., 0.8, &
2603 (/12.,0./), (/34.,2./)),
"Salt stratified Linearized Alpha/Beta"))
2604 if (.not. ndiff_unit_tests_discontinuous)
write(stdout,*)
'Pass'
2606 end function ndiff_unit_tests_discontinuous
2609 logical function test_fv_diff(verbose, hkm1, hk, hkp1, Skm1, Sk, Skp1, Ptrue, title)
2610 logical,
intent(in) :: verbose
2611 real,
intent(in) :: hkm1
2612 real,
intent(in) :: hk
2613 real,
intent(in) :: hkp1
2614 real,
intent(in) :: skm1
2615 real,
intent(in) :: sk
2616 real,
intent(in) :: skp1
2617 real,
intent(in) :: ptrue
2618 character(len=*),
intent(in) :: title
2624 pret = fv_diff(hkm1, hk, hkp1, skm1, sk, skp1)
2625 test_fv_diff = (pret /= ptrue)
2627 if (test_fv_diff .or. verbose)
then
2629 if (test_fv_diff) stdunit = stderr
2630 write(stdunit,
'(a)') title
2631 if (test_fv_diff)
then
2632 write(stdunit,
'(2(x,a,f20.16),x,a)')
'pRet=',pret,
'pTrue=',ptrue,
'WRONG!'
2634 write(stdunit,
'(2(x,a,f20.16))')
'pRet=',pret,
'pTrue=',ptrue
2638 end function test_fv_diff
2641 logical function test_fvlsq_slope(verbose, hkm1, hk, hkp1, Skm1, Sk, Skp1, Ptrue, title)
2642 logical,
intent(in) :: verbose
2643 real,
intent(in) :: hkm1
2644 real,
intent(in) :: hk
2645 real,
intent(in) :: hkp1
2646 real,
intent(in) :: skm1
2647 real,
intent(in) :: sk
2648 real,
intent(in) :: skp1
2649 real,
intent(in) :: ptrue
2650 character(len=*),
intent(in) :: title
2656 pret = fvlsq_slope(hkm1, hk, hkp1, skm1, sk, skp1)
2657 test_fvlsq_slope = (pret /= ptrue)
2659 if (test_fvlsq_slope .or. verbose)
then
2661 if (test_fvlsq_slope) stdunit = stderr
2662 write(stdunit,
'(a)') title
2663 if (test_fvlsq_slope)
then
2664 write(stdunit,
'(2(x,a,f20.16),x,a)')
'pRet=',pret,
'pTrue=',ptrue,
'WRONG!'
2666 write(stdunit,
'(2(x,a,f20.16))')
'pRet=',pret,
'pTrue=',ptrue
2670 end function test_fvlsq_slope
2673 logical function test_ifndp(verbose, rhoNeg, Pneg, rhoPos, Ppos, Ptrue, title)
2674 logical,
intent(in) :: verbose
2675 real,
intent(in) :: rhoneg
2676 real,
intent(in) :: pneg
2677 real,
intent(in) :: rhopos
2678 real,
intent(in) :: ppos
2679 real,
intent(in) :: ptrue
2680 character(len=*),
intent(in) :: title
2686 pret = interpolate_for_nondim_position(rhoneg, pneg, rhopos, ppos)
2687 test_ifndp = (pret /= ptrue)
2689 if (test_ifndp .or. verbose)
then
2691 if (test_ifndp) stdunit = stderr
2692 write(stdunit,
'(a)') title
2693 if (test_ifndp)
then
2694 write(stdunit,
'(4(x,a,f20.16),2(x,a,1pe22.15),x,a)') &
2695 'r1=',rhoneg,
'p1=',pneg,
'r2=',rhopos,
'p2=',ppos,
'pRet=',pret,
'pTrue=',ptrue,
'WRONG!'
2697 write(stdunit,
'(4(x,a,f20.16),2(x,a,1pe22.15))') &
2698 'r1=',rhoneg,
'p1=',pneg,
'r2=',rhopos,
'p2=',ppos,
'pRet=',pret,
'pTrue=',ptrue
2702 end function test_ifndp
2705 logical function test_data1d(verbose, nk, Po, Ptrue, title)
2706 logical,
intent(in) :: verbose
2707 integer,
intent(in) :: nk
2708 real,
dimension(nk),
intent(in) :: po
2709 real,
dimension(nk),
intent(in) :: ptrue
2710 character(len=*),
intent(in) :: title
2713 integer :: k, stdunit
2715 test_data1d = .false.
2717 if (po(k) /= ptrue(k)) test_data1d = .true.
2720 if (test_data1d .or. verbose)
then
2722 if (test_data1d) stdunit = stderr
2723 write(stdunit,
'(a)') title
2725 if (po(k) /= ptrue(k))
then
2726 test_data1d = .true.
2727 write(stdunit,
'(a,i2,2(x,a,f20.16),x,a,1pe22.15,x,a)') &
2728 'k=',k,
'Po=',po(k),
'Ptrue=',ptrue(k),
'err=',po(k)-ptrue(k),
'WRONG!'
2731 write(stdunit,
'(a,i2,2(x,a,f20.16),x,a,1pe22.15)') &
2732 'k=',k,
'Po=',po(k),
'Ptrue=',ptrue(k),
'err=',po(k)-ptrue(k)
2737 end function test_data1d
2740 logical function test_data1di(verbose, nk, Po, Ptrue, title)
2741 logical,
intent(in) :: verbose
2742 integer,
intent(in) :: nk
2743 integer,
dimension(nk),
intent(in) :: po
2744 integer,
dimension(nk),
intent(in) :: ptrue
2745 character(len=*),
intent(in) :: title
2748 integer :: k, stdunit
2750 test_data1di = .false.
2752 if (po(k) /= ptrue(k)) test_data1di = .true.
2755 if (test_data1di .or. verbose)
then
2757 if (test_data1di) stdunit = stderr
2758 write(stdunit,
'(a)') title
2760 if (po(k) /= ptrue(k))
then
2761 test_data1di = .true.
2762 write(stdunit,
'(a,i2,2(x,a,i5),x,a)')
'k=',k,
'Io=',po(k),
'Itrue=',ptrue(k),
'WRONG!'
2765 write(stdunit,
'(a,i2,2(x,a,i5))')
'k=',k,
'Io=',po(k),
'Itrue=',ptrue(k)
2770 end function test_data1di
2774 logical function test_nsp(verbose, ns, KoL, KoR, pL, pR, hEff, KoL0, KoR0, pL0, pR0, hEff0, title)
2775 logical,
intent(in) :: verbose
2776 integer,
intent(in) :: ns
2777 integer,
dimension(ns),
intent(in) :: kol
2778 integer,
dimension(ns),
intent(in) :: kor
2779 real,
dimension(ns),
intent(in) :: pl
2780 real,
dimension(ns),
intent(in) :: pr
2781 real,
dimension(ns-1),
intent(in) :: heff
2782 integer,
dimension(ns),
intent(in) :: kol0
2783 integer,
dimension(ns),
intent(in) :: kor0
2784 real,
dimension(ns),
intent(in) :: pl0
2785 real,
dimension(ns),
intent(in) :: pr0
2786 real,
dimension(ns-1),
intent(in) :: heff0
2787 character(len=*),
intent(in) :: title
2790 integer :: k, stdunit
2791 logical :: this_row_failed
2795 test_nsp = test_nsp .or. compare_nsp_row(kol(k), kor(k), pl(k), pr(k), kol0(k), kor0(k), pl0(k), pr0(k))
2797 if (heff(k) /= heff0(k)) test_nsp = .true.
2801 if (test_nsp .or. verbose)
then
2803 if (test_nsp) stdunit = stderr
2804 write(stdunit,
'(a)') title
2806 this_row_failed = compare_nsp_row(kol(k), kor(k), pl(k), pr(k), kol0(k), kor0(k), pl0(k), pr0(k))
2807 if (this_row_failed)
then
2808 write(stdunit,10) k,kol(k),pl(k),kor(k),pr(k),
' <-- WRONG!'
2809 write(stdunit,10) k,kol0(k),pl0(k),kor0(k),pr0(k),
' <-- should be this'
2811 write(stdunit,10) k,kol(k),pl(k),kor(k),pr(k)
2814 if (heff(k) /= heff0(k))
then
2815 write(stdunit,
'(i3,8x,"layer hEff =",2(f20.16,a))') k,heff(k),
" .neq. ",heff0(k),
' <-- WRONG!'
2817 write(stdunit,
'(i3,8x,"layer hEff =",f20.16)') k,heff(k)
2822 if (test_nsp)
call mom_error(fatal,
"test_nsp failed")
2824 10
format(
"ks=",i3,
" kL=",i3,
" pL=",f20.16,
" kR=",i3,
" pR=",f20.16,a)
2825 end function test_nsp
2828 logical function compare_nsp_row(KoL, KoR, pL, pR, KoL0, KoR0, pL0, pR0)
2829 integer,
intent(in) :: kol
2830 integer,
intent(in) :: kor
2831 real,
intent(in) :: pl
2832 real,
intent(in) :: pr
2833 integer,
intent(in) :: kol0
2834 integer,
intent(in) :: kor0
2835 real,
intent(in) :: pl0
2836 real,
intent(in) :: pr0
2838 compare_nsp_row = .false.
2839 if (kol /= kol0) compare_nsp_row = .true.
2840 if (kor /= kor0) compare_nsp_row = .true.
2841 if (pl /= pl0) compare_nsp_row = .true.
2842 if (pr /= pr0) compare_nsp_row = .true.
2843 end function compare_nsp_row
2846 logical function test_rnp(expected_pos, test_pos, title)
2847 real,
intent(in) :: expected_pos
2848 real,
intent(in) :: test_pos
2849 character(len=*),
intent(in) :: title
2854 test_rnp = abs(expected_pos - test_pos) > 2*epsilon(test_pos)
2856 write(stdunit,
'(A, f20.16, " .neq. ", f20.16, " <-- WRONG")') title, expected_pos, test_pos
2858 write(stdunit,
'(A, f20.16, " == ", f20.16)') title, expected_pos, test_pos
2860 end function test_rnp
2862 subroutine neutral_diffusion_end(CS)
2865 if (
associated(cs))
deallocate(cs)
2867 end subroutine neutral_diffusion_end