MOM6
coord_slight.F90
1 !> Regrid columns for the SLight coordinate
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_error_handler, only : mom_error, fatal
9 use regrid_interp, only : interp_cs_type, regridding_set_ppolys
10 use regrid_interp, only : nr_iterations, nr_tolerance, degree_max
11 
12 implicit none ; private
13 
14 !> Control structure containing required parameters for the SLight coordinate
15 type, public :: slight_cs ; private
16 
17  !> Number of layers/levels
18  integer :: nk
19 
20  !> Minimum thickness allowed when building the new grid through regridding [H ~> m or kg m-2]
21  real :: min_thickness
22 
23  !> Reference pressure for potential density calculations [R L2 T-2 ~> Pa]
24  real :: ref_pressure
25 
26  !> Fraction (between 0 and 1) of compressibility to add to potential density
27  !! profiles when interpolating for target grid positions. [nondim]
28  real :: compressibility_fraction
29 
30  ! The following 4 parameters were introduced for use with the SLight coordinate:
31  !> Depth over which to average to determine the mixed layer potential density [H ~> m or kg m-2]
32  real :: rho_ml_avg_depth
33 
34  !> Number of layers to offset the mixed layer density to find resolved stratification [nondim]
35  real :: nlay_ml_offset
36 
37  !> The number of fixed-thickness layers at the top of the model
38  integer :: nz_fixed_surface = 2
39 
40  !> The fixed resolution in the topmost SLight_nkml_min layers [H ~> m or kg m-2]
41  real :: dz_ml_min
42 
43  !> If true, detect regions with much weaker stratification in the coordinate
44  !! than based on in-situ density, and use a stretched coordinate there.
45  logical :: fix_haloclines = .false.
46 
47  !> A length scale over which to filter T & S when looking for spuriously
48  !! unstable water mass profiles [H ~> m or kg m-2].
49  real :: halocline_filter_length
50 
51  !> A value of the stratification ratio that defines a problematic halocline region [nondim].
52  real :: halocline_strat_tol
53 
54  !> Nominal density of interfaces [R ~> kg m-3].
55  real, allocatable, dimension(:) :: target_density
56 
57  !> Maximum depths of interfaces [H ~> m or kg m-2].
58  real, allocatable, dimension(:) :: max_interface_depths
59 
60  !> Maximum thicknesses of layers [H ~> m or kg m-2].
61  real, allocatable, dimension(:) :: max_layer_thickness
62 
63  !> Interpolation control structure
64  type(interp_cs_type) :: interp_cs
65 end type slight_cs
66 
67 public init_coord_slight, set_slight_params, build_slight_column, end_coord_slight
68 
69 contains
70 
71 !> Initialise a slight_CS with pointers to parameters
72 subroutine init_coord_slight(CS, nk, ref_pressure, target_density, interp_CS, m_to_H)
73  type(slight_cs), pointer :: cs !< Unassociated pointer to hold the control structure
74  integer, intent(in) :: nk !< Number of layers in the grid
75  real, intent(in) :: ref_pressure !< Coordinate reference pressure [R L2 T-2 ~> Pa]
76  real, dimension(:), intent(in) :: target_density !< Nominal density of interfaces [R ~> kg m-3]
77  type(interp_cs_type), intent(in) :: interp_cs !< Controls for interpolation
78  real, optional, intent(in) :: m_to_h !< A conversion factor from m to the units of thicknesses
79 
80  real :: m_to_h_rescale ! A unit conversion factor.
81 
82  if (associated(cs)) call mom_error(fatal, "init_coord_slight: CS already associated!")
83  allocate(cs)
84  allocate(cs%target_density(nk+1))
85 
86  m_to_h_rescale = 1.0 ; if (present(m_to_h)) m_to_h_rescale = m_to_h
87 
88  cs%nk = nk
89  cs%ref_pressure = ref_pressure
90  cs%target_density(:) = target_density(:)
91  cs%interp_CS = interp_cs
92 
93  ! Set real parameter default values
94  cs%compressibility_fraction = 0. ! Nondim.
95  cs%Rho_ML_avg_depth = 1.0 * m_to_h_rescale
96  cs%nlay_ml_offset = 2.0 ! Nondim.
97  cs%dz_ml_min = 1.0 * m_to_h_rescale
98  cs%halocline_filter_length = 2.0 * m_to_h_rescale
99  cs%halocline_strat_tol = 0.25 ! Nondim.
100 
101 end subroutine init_coord_slight
102 
103 !> This subroutine deallocates memory in the control structure for the coord_slight module
104 subroutine end_coord_slight(CS)
105  type(slight_cs), pointer :: cs !< Coordinate control structure
106 
107  ! nothing to do
108  if (.not. associated(cs)) return
109  deallocate(cs%target_density)
110  deallocate(cs)
111 end subroutine end_coord_slight
112 
113 !> This subroutine can be used to set the parameters for the coord_slight module
114 subroutine set_slight_params(CS, max_interface_depths, max_layer_thickness, &
115  min_thickness, compressibility_fraction, dz_ml_min, &
116  nz_fixed_surface, Rho_ML_avg_depth, nlay_ML_offset, fix_haloclines, &
117  halocline_filter_length, halocline_strat_tol, interp_CS)
118  type(slight_cs), pointer :: cs !< Coordinate control structure
119  real, dimension(:), &
120  optional, intent(in) :: max_interface_depths !< Maximum depths of interfaces [H ~> m or kg m-2]
121  real, dimension(:), &
122  optional, intent(in) :: max_layer_thickness !< Maximum thicknesses of layers [H ~> m or kg m-2]
123  real, optional, intent(in) :: min_thickness !< Minimum thickness allowed when building the
124  !! new grid through regridding [H ~> m or kg m-2]
125  real, optional, intent(in) :: compressibility_fraction !< Fraction (between 0 and 1) of
126  !! compressibility to add to potential density profiles when
127  !! interpolating for target grid positions. [nondim]
128  real, optional, intent(in) :: dz_ml_min !< The fixed resolution in the topmost
129  !! SLight_nkml_min layers [H ~> m or kg m-2]
130  integer, optional, intent(in) :: nz_fixed_surface !< The number of fixed-thickness layers at the
131  !! top of the model
132  real, optional, intent(in) :: rho_ml_avg_depth !< Depth over which to average to determine
133  !! the mixed layer potential density [H ~> m or kg m-2]
134  real, optional, intent(in) :: nlay_ml_offset !< Number of layers to offset the mixed layer
135  !! density to find resolved stratification [nondim]
136  logical, optional, intent(in) :: fix_haloclines !< If true, detect regions with much weaker than
137  !! based on in-situ density, and use a stretched coordinate there.
138  real, optional, intent(in) :: halocline_filter_length !< A length scale over which to filter T & S
139  !! when looking for spuriously unstable water mass profiles [H ~> m or kg m-2].
140  real, optional, intent(in) :: halocline_strat_tol !< A value of the stratification ratio that
141  !! defines a problematic halocline region [nondim].
142  type(interp_cs_type), &
143  optional, intent(in) :: interp_cs !< Controls for interpolation
144 
145  if (.not. associated(cs)) call mom_error(fatal, "set_slight_params: CS not associated")
146 
147  if (present(max_interface_depths)) then
148  if (size(max_interface_depths) /= cs%nk+1) &
149  call mom_error(fatal, "set_slight_params: max_interface_depths inconsistent size")
150  allocate(cs%max_interface_depths(cs%nk+1))
151  cs%max_interface_depths(:) = max_interface_depths(:)
152  endif
153 
154  if (present(max_layer_thickness)) then
155  if (size(max_layer_thickness) /= cs%nk) &
156  call mom_error(fatal, "set_slight_params: max_layer_thickness inconsistent size")
157  allocate(cs%max_layer_thickness(cs%nk))
158  cs%max_layer_thickness(:) = max_layer_thickness(:)
159  endif
160 
161  if (present(min_thickness)) cs%min_thickness = min_thickness
162  if (present(compressibility_fraction)) cs%compressibility_fraction = compressibility_fraction
163 
164  if (present(dz_ml_min)) cs%dz_ml_min = dz_ml_min
165  if (present(nz_fixed_surface)) cs%nz_fixed_surface = nz_fixed_surface
166  if (present(rho_ml_avg_depth)) cs%Rho_ML_avg_depth = rho_ml_avg_depth
167  if (present(nlay_ml_offset)) cs%nlay_ML_offset = nlay_ml_offset
168  if (present(fix_haloclines)) cs%fix_haloclines = fix_haloclines
169  if (present(halocline_filter_length)) cs%halocline_filter_length = halocline_filter_length
170  if (present(halocline_strat_tol)) then
171  if (halocline_strat_tol > 1.0) call mom_error(fatal, "set_slight_params: "//&
172  "HALOCLINE_STRAT_TOL must not exceed 1.0.")
173  cs%halocline_strat_tol = halocline_strat_tol
174  endif
175 
176  if (present(interp_cs)) cs%interp_CS = interp_cs
177 end subroutine set_slight_params
178 
179 !> Build a SLight coordinate column
180 subroutine build_slight_column(CS, eqn_of_state, H_to_pres, H_subroundoff, &
181  nz, depth, h_col, T_col, S_col, p_col, z_col, z_col_new, &
182  h_neglect, h_neglect_edge)
183  type(slight_cs), intent(in) :: cs !< Coordinate control structure
184  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
185  real, intent(in) :: h_to_pres !< A conversion factor from thicknesses to
186  !! scaled pressure [R L2 T-2 H-1 ~> Pa m-1 or Pa m2 kg-1]
187  real, intent(in) :: h_subroundoff !< GV%H_subroundoff
188  integer, intent(in) :: nz !< Number of levels
189  real, intent(in) :: depth !< Depth of ocean bottom (positive [H ~> m or kg m-2])
190  real, dimension(nz), intent(in) :: t_col !< T for column
191  real, dimension(nz), intent(in) :: s_col !< S for column
192  real, dimension(nz), intent(in) :: h_col !< Layer thicknesses [H ~> m or kg m-2]
193  real, dimension(nz), intent(in) :: p_col !< Layer center pressure [R L2 T-2 ~> Pa]
194  real, dimension(nz+1), intent(in) :: z_col !< Interface positions relative to the surface [H ~> m or kg m-2]
195  real, dimension(nz+1), intent(inout) :: z_col_new !< Absolute positions of interfaces [H ~> m or kg m-2]
196  real, optional, intent(in) :: h_neglect !< A negligibly small width for the purpose of
197  !! cell reconstructions [H ~> m or kg m-2].
198  real, optional, intent(in) :: h_neglect_edge !< A negligibly small width for the purpose
199  !! of edge value calculations [H ~> m or kg m-2].
200  ! Local variables
201  real, dimension(nz) :: rho_col ! Layer densities [R ~> kg m-3]
202  real, dimension(nz) :: t_f, s_f ! Filtered layer temperature [degC] and salinity [ppt]
203  logical, dimension(nz+1) :: reliable ! If true, this interface is in a reliable position.
204  real, dimension(nz+1) :: t_int, s_int ! Temperature [degC] and salinity [ppt] interpolated to interfaces.
205  real, dimension(nz+1) :: rho_tmp ! A temporary density [R ~> kg m-3]
206  real, dimension(nz+1) :: drho_dp ! The partial derivative of density with pressure [T2 L-2 ~> kg m-3 Pa-1]
207  real, dimension(nz+1) :: p_is, p_r ! Pressures [R L2 T-2 ~> Pa]
208  real, dimension(nz+1) :: drhois_dt ! The partial derivative of in situ density with temperature
209  ! in [R degC-1 ~> kg m-3 degC-1]
210  real, dimension(nz+1) :: drhois_ds ! The partial derivative of in situ density with salinity
211  ! in [R ppt-1 ~> kg m-3 ppt-1]
212  real, dimension(nz+1) :: drhor_dt ! The partial derivative of reference density with temperature
213  ! in [R degC-1 ~> kg m-3 degC-1]
214  real, dimension(nz+1) :: drhor_ds ! The partial derivative of reference density with salinity
215  ! in [R ppt-1 ~> kg m-3 ppt-1]
216  real, dimension(nz+1) :: strat_rat
217  real :: h_to_cpa ! A conversion factor from thicknesses to the compressibility fraction times
218  ! the units of pressure [R L2 T-2 H-1 ~> Pa m-1 or Pa m2 kg-1]
219  real :: dris, drr ! In situ and reference density differences [R ~> kg m-3]
220  real :: fn_now, i_hstol, fn_zero_val ! Nondimensional variables [nondim]
221  real :: z_int_unst ! The depth where the stratification allows the interior grid to start [H ~> m or kg m-2]
222  real :: dz ! A uniform layer thickness in very shallow water [H ~> m or kg m-2].
223  real :: dz_ur ! The total thickness of an unstable region [H ~> m or kg m-2].
224  real :: wgt, cowgt ! A weight and its complement [nondim].
225  real :: rho_ml_av ! The average potential density in a near-surface region [R ~> kg m-3].
226  real :: h_ml_av ! A thickness to try to use in taking the near-surface average [H ~> m or kg m-2].
227  real :: rho_x_z ! A cumulative integral of a density [R H ~> kg m-2 or kg2 m-5].
228  real :: z_wt ! The thickness actually used in taking the near-surface average [H ~> m or kg m-2].
229  real :: k_interior ! The (real) value of k where the interior grid starts [nondim].
230  real :: k_int2 ! The (real) value of k where the interior grid starts [nondim].
231  real :: z_interior ! The depth where the interior grid starts [H ~> m or kg m-2].
232  real :: z_ml_fix ! The depth at which the fixed-thickness near-surface layers end [H ~> m or kg m-2].
233  real :: dz_dk ! The thickness of layers between the fixed-thickness
234  ! near-surface layars and the interior [H ~> m or kg m-2].
235  real :: lfilt ! A filtering lengthscale [H ~> m or kg m-2].
236  logical :: maximum_depths_set ! If true, the maximum depths of interface have been set.
237  logical :: maximum_h_set ! If true, the maximum layer thicknesses have been set.
238  real :: k2_used, k2here, dz_sum, z_max
239  integer :: k2
240  real :: h_tr, b_denom_1, b1, d1 ! Temporary variables used by the tridiagonal solver.
241  real, dimension(nz) :: c1 ! Temporary variables used by the tridiagonal solver.
242  integer :: kur1, kur2 ! The indicies at the top and bottom of an unreliable region.
243  integer :: kur_ss ! The index to start with in the search for the next unstable region.
244  integer :: i, j, k, nkml
245 
246  maximum_depths_set = allocated(cs%max_interface_depths)
247  maximum_h_set = allocated(cs%max_layer_thickness)
248 
249  if (z_col(nz+1) - z_col(1) < nz*cs%min_thickness) then
250  ! This is a nearly massless total depth, so distribute the water evenly.
251  dz = (z_col(nz+1) - z_col(1)) / real(nz)
252  do k=2,nz ; z_col_new(k) = z_col(1) + dz*real(k-1) ; enddo
253  else
254  call calculate_density(t_col, s_col, p_col, rho_col, eqn_of_state)
255 
256  ! Find the locations of the target potential densities, flagging
257  ! locations in apparently unstable regions as not reliable.
258  call rho_interfaces_col(rho_col, h_col, z_col, cs%target_density, nz, &
259  z_col_new, cs, reliable, debug=.true., &
260  h_neglect=h_neglect, h_neglect_edge=h_neglect_edge)
261 
262  ! Ensure that the interfaces are at least CS%min_thickness apart.
263  if (cs%min_thickness > 0.0) then
264  ! Move down interfaces below overly thin layers.
265  do k=2,nz ; if (z_col_new(k) < z_col_new(k-1) + cs%min_thickness) then
266  z_col_new(k) = z_col_new(k-1) + cs%min_thickness
267  endif ; enddo
268  ! Now move up any interfaces that are too close to the bottom.
269  do k=nz,2,-1 ; if (z_col_new(k) > z_col_new(k+1) - cs%min_thickness) then
270  z_col_new(k) = z_col_new(k+1) - cs%min_thickness
271  else
272  exit ! No more interfaces can be too close to the bottom.
273  endif ; enddo
274  endif
275 
276  ! Fix up the unreliable regions.
277  kur_ss = 2 ! reliable(1) and reliable(nz+1) must always be true.
278  do
279  ! Search for the uppermost unreliable interface postion.
280  kur1 = nz+2
281  do k=kur_ss,nz ; if (.not.reliable(k)) then
282  kur1 = k ; exit
283  endif ; enddo
284  if (kur1 > nz) exit ! Everything is now reliable.
285 
286  kur2 = kur1-1 ! For error checking.
287  do k=kur1+1,nz+1 ; if (reliable(k)) then
288  kur2 = k-1 ; kur_ss = k ; exit
289  endif ; enddo
290  if (kur2 < kur1) call mom_error(fatal, "Bad unreliable range.")
291 
292  dz_ur = z_col_new(kur2+1) - z_col_new(kur1-1)
293  ! drho = CS%target_density(kur2+1) - CS%target_density(kur1-1)
294  ! Perhaps reset the wgt and cowgt depending on how bad the old interface
295  ! locations were.
296  wgt = 1.0 ; cowgt = 0.0 ! = 1.0-wgt
297  do k=kur1,kur2
298  z_col_new(k) = cowgt*z_col_new(k) + &
299  wgt * (z_col_new(kur1-1) + dz_ur*(k - (kur1-1)) / ((kur2 - kur1) + 2))
300  enddo
301  enddo
302 
303  ! Determine which interfaces are in the s-space region and the depth extent
304  ! of this region.
305  z_wt = 0.0 ; rho_x_z = 0.0
306  h_ml_av = cs%Rho_ml_avg_depth
307  do k=1,nz
308  if (z_wt + h_col(k) >= h_ml_av) then
309  rho_x_z = rho_x_z + rho_col(k) * (h_ml_av - z_wt)
310  z_wt = h_ml_av
311  exit
312  else
313  rho_x_z = rho_x_z + rho_col(k) * h_col(k)
314  z_wt = z_wt + h_col(k)
315  endif
316  enddo
317  if (z_wt > 0.0) rho_ml_av = rho_x_z / z_wt
318 
319  nkml = cs%nz_fixed_surface
320  ! Find the interface that matches rho_ml_av.
321  if (rho_ml_av <= cs%target_density(nkml)) then
322  k_interior = cs%nlay_ml_offset + real(nkml)
323  elseif (rho_ml_av > cs%target_density(nz+1)) then
324  k_interior = real(nz+1)
325  else ; do k=nkml,nz
326  if ((rho_ml_av >= cs%target_density(k)) .and. &
327  (rho_ml_av < cs%target_density(k+1))) then
328  k_interior = (cs%nlay_ml_offset + k) + &
329  (rho_ml_av - cs%target_density(k)) / &
330  (cs%target_density(k+1) - cs%target_density(k))
331  exit
332  endif
333  enddo ; endif
334  if (k_interior > real(nz+1)) k_interior = real(nz+1)
335 
336  ! Linearly interpolate to find z_interior. This could be made more sophisticated.
337  k = int(ceiling(k_interior))
338  z_interior = (k-k_interior)*z_col_new(k-1) + (1.0+(k_interior-k))*z_col_new(k)
339 
340  if (cs%fix_haloclines) then
341  ! ! Identify regions above the reference pressure where the chosen
342  ! ! potential density significantly underestimates the actual
343  ! ! stratification, and use these to find a second estimate of
344  ! ! z_int_unst and k_interior.
345 
346  if (cs%halocline_filter_length > 0.0) then
347  lfilt = cs%halocline_filter_length
348 
349  ! Filter the temperature and salnity with a fixed lengthscale.
350  h_tr = h_col(1) + h_subroundoff
351  b1 = 1.0 / (h_tr + lfilt) ; d1 = h_tr * b1
352  t_f(1) = (b1*h_tr)*t_col(1) ; s_f(1) = (b1*h_tr)*s_col(1)
353  do k=2,nz
354  c1(k) = lfilt * b1
355  h_tr = h_col(k) + h_subroundoff ; b_denom_1 = h_tr + d1*lfilt
356  b1 = 1.0 / (b_denom_1 + lfilt) ; d1 = b_denom_1 * b1
357  t_f(k) = b1 * (h_tr*t_col(k) + lfilt*t_f(k-1))
358  s_f(k) = b1 * (h_tr*s_col(k) + lfilt*s_f(k-1))
359  enddo
360  do k=nz-1,1,-1
361  t_f(k) = t_f(k) + c1(k+1)*t_f(k+1) ; s_f(k) = s_f(k) + c1(k+1)*s_f(k+1)
362  enddo
363  else
364  do k=1,nz ; t_f(k) = t_col(k) ; s_f(k) = s_col(k) ; enddo
365  endif
366 
367  t_int(1) = t_f(1) ; s_int(1) = s_f(1)
368  do k=2,nz
369  t_int(k) = 0.5*(t_f(k-1) + t_f(k)) ; s_int(k) = 0.5*(s_f(k-1) + s_f(k))
370  p_is(k) = z_col(k) * h_to_pres
371  p_r(k) = cs%ref_pressure + cs%compressibility_fraction * ( p_is(k) - cs%ref_pressure )
372  enddo
373  t_int(nz+1) = t_f(nz) ; s_int(nz+1) = s_f(nz)
374  p_is(nz+1) = z_col(nz+1) * h_to_pres
375  call calculate_density_derivs(t_int, s_int, p_is, drhois_dt, drhois_ds, &
376  eqn_of_state, (/2,nz/) )
377  call calculate_density_derivs(t_int, s_int, p_r, drhor_dt, drhor_ds, &
378  eqn_of_state, (/2,nz/) )
379  if (cs%compressibility_fraction > 0.0) then
380  call calculate_compress(t_int, s_int, p_r(:), rho_tmp, drho_dp, 2, nz-1, eqn_of_state)
381  else
382  do k=2,nz ; drho_dp(k) = 0.0 ; enddo
383  endif
384 
385  h_to_cpa = cs%compressibility_fraction * h_to_pres
386  strat_rat(1) = 1.0
387  do k=2,nz
388  dris = drhois_dt(k) * (t_f(k) - t_f(k-1)) + &
389  drhois_ds(k) * (s_f(k) - s_f(k-1))
390  drr = (drhor_dt(k) * (t_f(k) - t_f(k-1)) + &
391  drhor_ds(k) * (s_f(k) - s_f(k-1))) + &
392  drho_dp(k) * (h_to_cpa*0.5*(h_col(k) + h_col(k-1)))
393 
394  if (dris <= 0.0) then
395  strat_rat(k) = 2.0 ! Maybe do this? => ; if (drR < 0.0) strat_rat(K) = -2.0
396  else
397  strat_rat(k) = 2.0*max(drr,0.0) / (dris + abs(drr))
398  endif
399  enddo
400  strat_rat(nz+1) = 1.0
401 
402  z_int_unst = 0.0 ; fn_now = 0.0
403  fn_zero_val = min(2.0*cs%halocline_strat_tol, &
404  0.5*(1.0 + cs%halocline_strat_tol))
405  if (cs%halocline_strat_tol > 0.0) then
406  ! Use Adcroft's reciprocal rule.
407  i_hstol = 0.0 ; if (fn_zero_val - cs%halocline_strat_tol > 0.0) &
408  i_hstol = 1.0 / (fn_zero_val - cs%halocline_strat_tol)
409  do k=nz,1,-1 ; if (cs%ref_pressure > p_is(k+1)) then
410  z_int_unst = z_int_unst + fn_now * h_col(k)
411  if (strat_rat(k) <= fn_zero_val) then
412  if (strat_rat(k) <= cs%halocline_strat_tol) then ; fn_now = 1.0
413  else
414  fn_now = max(fn_now, (fn_zero_val - strat_rat(k)) * i_hstol)
415  endif
416  endif
417  endif ; enddo
418  else
419  do k=nz,1,-1 ; if (cs%ref_pressure > p_is(k+1)) then
420  z_int_unst = z_int_unst + fn_now * h_col(k)
421  if (strat_rat(k) <= cs%halocline_strat_tol) fn_now = 1.0
422  endif ; enddo
423  endif
424 
425  if (z_interior < z_int_unst) then
426  ! Find a second estimate of the extent of the s-coordinate region.
427  kur1 = max(int(ceiling(k_interior)),2)
428  if (z_col_new(kur1-1) < z_interior) then
429  k_int2 = kur1
430  do k = kur1,nz+1 ; if (z_col_new(k) >= z_int_unst) then
431  ! This is linear interpolation again.
432  if (z_col_new(k-1) >= z_int_unst) &
433  call mom_error(fatal,"build_grid_SLight, bad halocline structure.")
434  k_int2 = real(k-1) + (z_int_unst - z_col_new(k-1)) / &
435  (z_col_new(k) - z_col_new(k-1))
436  exit
437  endif ; enddo
438  if (z_col_new(nz+1) < z_int_unst) then
439  ! This should be unnecessary.
440  z_int_unst = z_col_new(nz+1) ; k_int2 = real(nz+1)
441  endif
442 
443  ! Now take the larger values.
444  if (k_int2 > k_interior) then
445  k_interior = k_int2 ; z_interior = z_int_unst
446  endif
447  endif
448  endif
449  endif ! fix_haloclines
450 
451  z_col_new(1) = 0.0
452  do k=2,nkml+1
453  z_col_new(k) = min((k-1)*cs%dz_ml_min, &
454  z_col_new(nz+1) - cs%min_thickness*(nz+1-k))
455  enddo
456  z_ml_fix = z_col_new(nkml+1)
457  if (z_interior > z_ml_fix) then
458  dz_dk = (z_interior - z_ml_fix) / (k_interior - (nkml+1))
459  do k=nkml+2,int(floor(k_interior))
460  z_col_new(k) = z_ml_fix + dz_dk * (k - (nkml+1))
461  enddo
462  else ! The fixed-thickness z-region penetrates into the interior.
463  do k=nkml+2,nz
464  if (z_col_new(k) <= z_col_new(cs%nz_fixed_surface+1)) then
465  z_col_new(k) = z_col_new(cs%nz_fixed_surface+1)
466  else ; exit ; endif
467  enddo
468  endif
469 
470  if (maximum_depths_set .and. maximum_h_set) then ; do k=2,nz
471  ! The loop bounds are 2 & nz so the top and bottom interfaces do not move.
472  ! Recall that z_col_new is positive downward.
473  z_col_new(k) = min(z_col_new(k), cs%max_interface_depths(k), &
474  z_col_new(k-1) + cs%max_layer_thickness(k-1))
475  enddo ; elseif (maximum_depths_set) then ; do k=2,nz
476  z_col_new(k) = min(z_col_new(k), cs%max_interface_depths(k))
477  enddo ; elseif (maximum_h_set) then ; do k=2,nz
478  z_col_new(k) = min(z_col_new(k), z_col_new(k-1) + cs%max_layer_thickness(k-1))
479  enddo ; endif
480 
481  endif ! Total thickness exceeds nz*CS%min_thickness.
482 
483 end subroutine build_slight_column
484 
485 !> Finds the new interface locations in a column of water that match the
486 !! prescribed target densities.
487 subroutine rho_interfaces_col(rho_col, h_col, z_col, rho_tgt, nz, z_col_new, &
488  CS, reliable, debug, h_neglect, h_neglect_edge)
489  integer, intent(in) :: nz !< Number of layers
490  real, dimension(nz), intent(in) :: rho_col !< Initial layer reference densities [R ~> kg m-3].
491  real, dimension(nz), intent(in) :: h_col !< Initial layer thicknesses [H ~> m or kg m-2].
492  real, dimension(nz+1), intent(in) :: z_col !< Initial interface heights [H ~> m or kg m-2].
493  real, dimension(nz+1), intent(in) :: rho_tgt !< Interface target densities.
494  real, dimension(nz+1), intent(inout) :: z_col_new !< New interface heights [H ~> m or kg m-2].
495  type(slight_cs), intent(in) :: CS !< Coordinate control structure
496  logical, dimension(nz+1), intent(inout) :: reliable !< If true, the interface positions
497  !! are well defined from a stable region.
498  logical, optional, intent(in) :: debug !< If present and true, do debugging checks.
499  real, optional, intent(in) :: h_neglect !< A negligibly small width for the purpose of
500  !! cell reconstructions [H ~> m or kg m-2]
501  real, optional, intent(in) :: h_neglect_edge !< A negligibly small width for the purpose
502  !! of edge value calculations [H ~> m or kg m-2]
503 
504  real, dimension(nz+1) :: ru_max_int ! The maximum and minimum densities in
505  real, dimension(nz+1) :: ru_min_int ! an unstable region around an interface [R ~> kg m-3].
506  real, dimension(nz) :: ru_max_lay ! The maximum and minimum densities in
507  real, dimension(nz) :: ru_min_lay ! an unstable region containing a layer [R ~> kg m-3].
508  real, dimension(nz,2) :: ppoly_i_E ! Edge value of polynomial [R ~> kg m-3]
509  real, dimension(nz,2) :: ppoly_i_S ! Edge slope of polynomial [R H-1 ~> kg m-4 or m-1]
510  real, dimension(nz,DEGREE_MAX+1) :: ppoly_i_coefficients ! Coefficients of polynomial [R ~> kg m-3]
511  logical, dimension(nz) :: unstable_lay ! If true, this layer is in an unstable region.
512  logical, dimension(nz+1) :: unstable_int ! If true, this interface is in an unstable region.
513  real :: rt ! The current target density [R ~> kg m-3].
514  real :: zf ! The fractional z-position within a layer of the target density [nondim].
515  real :: rfn ! The target density relative to the interpolated density [R ~> kg m-3]
516  real :: a(5) ! Coefficients of a local polynomial minus the target density [R ~> kg m-3].
517  real :: zf1, zf2 ! Two previous estimates of zf [nondim]
518  real :: rfn1, rfn2 ! Values of rfn at zf1 and zf2 [R ~> kg m-3]
519  real :: drfn_dzf ! The partial derivative of rfn with zf [R ~> kg m-3]
520  real :: sgn, delta_zf, zf_prev ! [nondim]
521  real :: tol ! The tolerance for convergence of zf [nondim]
522  logical :: k_found ! If true, the position has been found.
523  integer :: k_layer ! The index of the stable layer containing an interface.
524  integer :: ppoly_degree
525  integer :: k, k1, k1_min, itt, max_itt, m
526 
527  real :: z_sgn ! 1 or -1, depending on whether z increases with increasing K.
528  logical :: debugging
529 
530  debugging = .false. ; if (present(debug)) debugging = debug
531  max_itt = nr_iterations
532  tol = nr_tolerance
533 
534  z_sgn = 1.0 ; if ( z_col(1) > z_col(nz+1) ) z_sgn = -1.0
535  if (debugging) then
536  do k=1,nz
537  if (abs((z_col(k+1) - z_col(k)) - z_sgn*h_col(k)) > &
538  1.0e-14*(abs(z_col(k+1)) + abs(z_col(k)) + abs(h_col(k))) ) &
539  call mom_error(fatal, "rho_interfaces_col: Inconsistent z_col and h_col")
540  enddo
541  endif
542 
543  if ( z_col(1) == z_col(nz+1) ) then
544  ! This is a massless column!
545  do k=1,nz+1 ; z_col_new(k) = z_col(1) ; reliable(k) = .true. ; enddo
546  return
547  endif
548 
549  ! This sets up the piecewise polynomials based on the rho_col profile.
550  call regridding_set_ppolys(cs%interp_CS, rho_col, nz, h_col, ppoly_i_e, ppoly_i_s, &
551  ppoly_i_coefficients, ppoly_degree, h_neglect, h_neglect_edge)
552 
553  ! Determine the density ranges of unstably stratified segments.
554  ! Interfaces that start out in an unstably stratified segment can
555  ! only escape if they are outside of the bounds of that segment, and no
556  ! interfaces are ever mapped into an unstable segment.
557  unstable_int(1) = .false.
558  ru_max_int(1) = ppoly_i_e(1,1)
559 
560  unstable_lay(1) = (ppoly_i_e(1,1) > ppoly_i_e(1,2))
561  ru_max_lay(1) = max(ppoly_i_e(1,1), ppoly_i_e(1,2))
562 
563  do k=2,nz
564  unstable_int(k) = (ppoly_i_e(k-1,2) > ppoly_i_e(k,1))
565  ru_max_int(k) = max(ppoly_i_e(k-1,2), ppoly_i_e(k,1))
566  ru_min_int(k) = min(ppoly_i_e(k-1,2), ppoly_i_e(k,1))
567  if (unstable_int(k) .and. unstable_lay(k-1)) &
568  ru_max_int(k) = max(ru_max_lay(k-1), ru_max_int(k))
569 
570  unstable_lay(k) = (ppoly_i_e(k,1) > ppoly_i_e(k,2))
571  ru_max_lay(k) = max(ppoly_i_e(k,1), ppoly_i_e(k,2))
572  ru_min_lay(k) = min(ppoly_i_e(k,1), ppoly_i_e(k,2))
573  if (unstable_lay(k) .and. unstable_int(k)) &
574  ru_max_lay(k) = max(ru_max_int(k), ru_max_lay(k))
575  enddo
576  unstable_int(nz+1) = .false.
577  ru_min_int(nz+1) = ppoly_i_e(nz,2)
578 
579  do k=nz,1,-1
580  if (unstable_lay(k) .and. unstable_int(k+1)) &
581  ru_min_lay(k) = min(ru_min_int(k+1), ru_min_lay(k))
582 
583  if (unstable_int(k) .and. unstable_lay(k)) &
584  ru_min_int(k) = min(ru_min_lay(k), ru_min_int(k))
585  enddo
586 
587  z_col_new(1) = z_col(1) ; reliable(1) = .true.
588  k1_min = 1
589  do k=2,nz ! Find the locations of the various target densities for the interfaces.
590  rt = rho_tgt(k)
591  k_layer = -1
592  k_found = .false.
593 
594  ! Many light layers are found at the top, so start there.
595  if (rt <= ppoly_i_e(k1_min,1)) then
596  z_col_new(k) = z_col(k1_min)
597  k_found = .true.
598  ! Do not change k1_min for the next layer.
599  elseif (k1_min == nz+1) then
600  z_col_new(k) = z_col(nz+1)
601  else
602  ! Start with the previous location and search outward.
603  if (unstable_int(k) .and. (rt >= ru_min_int(k)) .and. (rt <= ru_max_int(k))) then
604  ! This interface started in an unstable region and should not move due to remapping.
605  z_col_new(k) = z_col(k) ; reliable(k) = .false.
606  k1_min = k ; k_found = .true.
607  elseif ((rt >= ppoly_i_e(k-1,2)) .and. (rt <= ppoly_i_e(k,1))) then
608  ! This interface is already in the right place and does not move.
609  z_col_new(k) = z_col(k) ; reliable(k) = .true.
610  k1_min = k ; k_found = .true.
611  elseif (rt < ppoly_i_e(k-1,2)) then ! Search upward
612  do k1=k-1,k1_min,-1
613  ! Check whether rt is in layer k.
614  if ((rt < ppoly_i_e(k1,2)) .and. (rt > ppoly_i_e(k1,1))) then
615  ! rt is in layer k.
616  k_layer = k1
617  k1_min = k1 ; k_found = .true. ; exit
618  elseif (unstable_lay(k1) .and. (rt >= ru_min_lay(k1)) .and. (rt <= ru_max_lay(k1))) then
619  ! rt would be found at unstable layer that it can not penetrate.
620  ! It is possible that this can never happen?
621  z_col_new(k) = z_col(k1+1) ; reliable(k) = .false.
622  k1_min = k1 ; k_found = .true. ; exit
623  endif
624  ! Check whether rt is at interface K.
625  if (k1 > 1) then ; if ((rt <= ppoly_i_e(k1,1)) .and. (rt >= ppoly_i_e(k1-1,2))) then
626  ! rt is at interface K1
627  z_col_new(k) = z_col(k1) ; reliable(k) = .true.
628  k1_min = k1 ; k_found = .true. ; exit
629  elseif (unstable_int(k1) .and. (rt >= ru_min_int(k1)) .and. (rt <= ru_max_int(k1))) then
630  ! rt would be found at an unstable interface that it can not pass.
631  ! It is possible that this can never happen?
632  z_col_new(k) = z_col(k1) ; reliable(k) = .false.
633  k1_min = k1 ; k_found = .true. ; exit
634  endif ; endif
635  enddo
636 
637  if (.not.k_found) then
638  ! This should not happen unless k1_min = 1.
639  if (k1_min < 2) then
640  z_col_new(k) = z_col(k1_min)
641  else
642  z_col_new(k) = z_col(k1_min)
643  endif
644  endif
645 
646  else ! Search downward
647  do k1=k,nz
648  if ((rt < ppoly_i_e(k1,2)) .and. (rt > ppoly_i_e(k1,1))) then
649  ! rt is in layer k.
650  k_layer = k1
651  k1_min = k1 ; k_found = .true. ; exit
652  elseif (unstable_lay(k1) .and. (rt >= ru_min_lay(k1)) .and. (rt <= ru_max_lay(k1))) then
653  ! rt would be found at unstable layer that it can not penetrate.
654  ! It is possible that this can never happen?
655  z_col_new(k) = z_col(k1)
656  reliable(k) = .false.
657  k1_min = k1 ; k_found = .true. ; exit
658  endif
659  if (k1 < nz) then ; if ((rt <= ppoly_i_e(k1+1,1)) .and. (rt >= ppoly_i_e(k1,2))) then
660  ! rt is at interface K1+1
661 
662  z_col_new(k) = z_col(k1+1) ; reliable(k) = .true.
663  k1_min = k1+1 ; k_found = .true. ; exit
664  elseif (unstable_int(k1+1) .and. (rt >= ru_min_int(k1+1)) .and. (rt <= ru_max_int(k1+1))) then
665  ! rt would be found at an unstable interface that it can not pass.
666  ! It is possible that this can never happen?
667  z_col_new(k) = z_col(k1+1)
668  reliable(k) = .false.
669  k1_min = k1+1 ; k_found = .true. ; exit
670  endif ; endif
671  enddo
672  if (.not.k_found) then
673  z_col_new(k) = z_col(nz+1)
674  if (rt >= ppoly_i_e(nz,2)) then
675  reliable(k) = .true.
676  else
677  reliable(k) = .false.
678  endif
679  endif
680  endif
681 
682  if (k_layer > 0) then ! The new location is inside of layer k_layer.
683  ! Note that this is coded assuming that this layer is stably stratified.
684  if (.not.(ppoly_i_e(k1,2) > ppoly_i_e(k1,1))) call mom_error(fatal, &
685  "build_grid_SLight: Erroneously searching for an interface in an unstratified layer.")
686 
687  ! Use the false position method to find the location (degree <= 1) or the first guess.
688  zf = (rt - ppoly_i_e(k1,1)) / (ppoly_i_e(k1,2) - ppoly_i_e(k1,1))
689 
690  if (ppoly_degree > 1) then ! Iterate to find the solution.
691  a(:) = 0.0 ; a(1) = ppoly_i_coefficients(k_layer,1) - rt
692  do m=2,ppoly_degree+1 ; a(m) = ppoly_i_coefficients(k_layer,m) ; enddo
693  ! Bracket the root.
694  zf1 = 0.0 ; rfn1 = a(1)
695  zf2 = 1.0 ; rfn2 = a(1) + (a(2) + (a(3) + (a(4) + a(5))))
696  if (rfn1 * rfn2 > 0.0) call mom_error(fatal, "build_grid_SLight: Bad bracketing.")
697 
698  do itt=1,max_itt
699  rfn = a(1) + zf*(a(2) + zf*(a(3) + zf*(a(4) + zf*a(5))))
700  ! Reset one of the ends of the bracket.
701  if (rfn * rfn1 > 0.0) then
702  zf1 = zf ; rfn1 = rfn
703  else
704  zf2 = zf ; rfn2 = rfn
705  endif
706  if (rfn1 == rfn2) exit
707 
708  drfn_dzf = (a(2) + zf*(2.0*a(3) + zf*(3.0*a(4) + zf*4.0*a(5))))
709  sgn = 1.0 ; if (drfn_dzf < 0.0) sgn = -1.0
710 
711  if ((sgn*(zf - rfn) >= zf1 * abs(drfn_dzf)) .and. &
712  (sgn*(zf - rfn) <= zf2 * abs(drfn_dzf))) then
713  delta_zf = -rfn / drfn_dzf
714  zf = zf + delta_zf
715  else ! Newton's method goes out of bounds, so use a false position method estimate
716  zf_prev = zf
717  zf = ( rfn2 * zf1 - rfn1 * zf2 ) / (rfn2 - rfn1)
718  delta_zf = zf - zf_prev
719  endif
720 
721  if (abs(delta_zf) < tol) exit
722  enddo
723  endif
724  z_col_new(k) = z_col(k_layer) + zf * z_sgn * h_col(k_layer)
725  reliable(k) = .true.
726  endif
727 
728  endif
729 
730  enddo
731  z_col_new(nz+1) = z_col(nz+1) ; reliable(nz+1) = .true.
732 
733 end subroutine rho_interfaces_col
734 
735 end module coord_slight
coord_slight::slight_cs
Control structure containing required parameters for the SLight coordinate.
Definition: coord_slight.F90:15
mom_eos::calculate_compress
Calculates the compressibility of water from T, S, and P.
Definition: MOM_EOS.F90:103
mom_eos
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
regrid_interp
Vertical interpolation for regridding.
Definition: regrid_interp.F90:2
mom_eos::eos_type
A control structure for the equation of state.
Definition: MOM_EOS.F90:108
mom_eos::calculate_density_derivs
Calculate the derivatives of density with temperature and salinity from T, S, and P.
Definition: MOM_EOS.F90:81
regrid_interp::interp_cs_type
Control structure for regrid_interp module.
Definition: regrid_interp.F90:23
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
coord_slight
Regrid columns for the SLight coordinate.
Definition: coord_slight.F90:2
mom_eos::calculate_density
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68