MOM6
coord_rho.F90
1 !> Regrid columns for the continuous isopycnal (rho) coordinate
2 module coord_rho
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_error_handler, only : mom_error, fatal
7 use mom_remapping, only : remapping_cs, remapping_core_h
9 use regrid_interp, only : interp_cs_type, build_and_interpolate_grid, degree_max
10 
11 implicit none ; private
12 
13 !> Control structure containing required parameters for the rho coordinate
14 type, public :: rho_cs ; private
15 
16  !> Number of layers
17  integer :: nk
18 
19  !> Minimum thickness allowed for layers, often in [H ~> m or kg m-2]
20  real :: min_thickness = 0.
21 
22  !> Reference pressure for density calculations [R L2 T-2 ~> Pa]
23  real :: ref_pressure
24 
25  !> If true, integrate for interface positions from the top downward.
26  !! If false, integrate from the bottom upward, as does the rest of the model.
27  logical :: integrate_downward_for_e = .false.
28 
29  !> Nominal density of interfaces [R ~> kg m-3]
30  real, allocatable, dimension(:) :: target_density
31 
32  !> Interpolation control structure
33  type(interp_cs_type) :: interp_cs
34 end type rho_cs
35 
36 public init_coord_rho, set_rho_params, build_rho_column, old_inflate_layers_1d, end_coord_rho
37 
38 contains
39 
40 !> Initialise a rho_CS with pointers to parameters
41 subroutine init_coord_rho(CS, nk, ref_pressure, target_density, interp_CS)
42  type(rho_cs), pointer :: CS !< Unassociated pointer to hold the control structure
43  integer, intent(in) :: nk !< Number of layers in the grid
44  real, intent(in) :: ref_pressure !< Coordinate reference pressure [R L2 T-2 ~> Pa]
45  real, dimension(:), intent(in) :: target_density !< Nominal density of interfaces [R ~> kg m-3]
46  type(interp_cs_type), intent(in) :: interp_CS !< Controls for interpolation
47 
48  if (associated(cs)) call mom_error(fatal, "init_coord_rho: CS already associated!")
49  allocate(cs)
50  allocate(cs%target_density(nk+1))
51 
52  cs%nk = nk
53  cs%ref_pressure = ref_pressure
54  cs%target_density(:) = target_density(:)
55  cs%interp_CS = interp_cs
56 
57 end subroutine init_coord_rho
58 
59 !> This subroutine deallocates memory in the control structure for the coord_rho module
60 subroutine end_coord_rho(CS)
61  type(rho_cs), pointer :: CS !< Coordinate control structure
62 
63  ! nothing to do
64  if (.not. associated(cs)) return
65  deallocate(cs%target_density)
66  deallocate(cs)
67 end subroutine end_coord_rho
68 
69 !> This subroutine can be used to set the parameters for the coord_rho module
70 subroutine set_rho_params(CS, min_thickness, integrate_downward_for_e, interp_CS)
71  type(rho_cs), pointer :: CS !< Coordinate control structure
72  real, optional, intent(in) :: min_thickness !< Minimum allowed thickness [H ~> m or kg m-2]
73  logical, optional, intent(in) :: integrate_downward_for_e !< If true, integrate for interface
74  !! positions from the top downward. If false, integrate
75  !! from the bottom upward, as does the rest of the model.
76 
77  type(interp_cs_type), optional, intent(in) :: interp_CS !< Controls for interpolation
78 
79  if (.not. associated(cs)) call mom_error(fatal, "set_rho_params: CS not associated")
80 
81  if (present(min_thickness)) cs%min_thickness = min_thickness
82  if (present(integrate_downward_for_e)) cs%integrate_downward_for_e = integrate_downward_for_e
83  if (present(interp_cs)) cs%interp_CS = interp_cs
84 end subroutine set_rho_params
85 
86 !> Build a rho coordinate column
87 !!
88 !! 1. Density profiles are calculated on the source grid.
89 !! 2. Positions of target densities (for interfaces) are found by interpolation.
90 subroutine build_rho_column(CS, nz, depth, h, T, S, eqn_of_state, z_interface, &
91  h_neglect, h_neglect_edge)
92  type(rho_cs), intent(in) :: CS !< coord_rho control structure
93  integer, intent(in) :: nz !< Number of levels on source grid (i.e. length of h, T, S)
94  real, intent(in) :: depth !< Depth of ocean bottom (positive downward) [H ~> m or kg m-2]
95  real, dimension(nz), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
96  real, dimension(nz), intent(in) :: T !< Temperature for source column [degC]
97  real, dimension(nz), intent(in) :: S !< Salinity for source column [ppt]
98  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
99  real, dimension(CS%nk+1), &
100  intent(inout) :: z_interface !< Absolute positions of interfaces
101  real, optional, intent(in) :: h_neglect !< A negligibly small width for the purpose
102  !! of cell reconstructions [H ~> m or kg m-2]
103  real, optional, intent(in) :: h_neglect_edge !< A negligibly small width for the purpose
104  !! of edge value calculations [H ~> m or kg m-2]
105 
106  ! Local variables
107  integer :: k, count_nonzero_layers
108  integer, dimension(nz) :: mapping
109  real, dimension(nz) :: pres ! Pressures used to calculate density [R L2 T-2 ~> Pa]
110  real, dimension(nz) :: h_nv ! Thicknesses of non-vanishing layers [H ~> m or kg m-2]
111  real, dimension(nz) :: densities ! Layer density [R ~> kg m-3]
112  real, dimension(nz+1) :: xTmp ! Temporary positions [H ~> m or kg m-2]
113  real, dimension(CS%nk) :: h_new ! New thicknesses [H ~> m or kg m-2]
114  real, dimension(CS%nk+1) :: x1 ! Interface heights [H ~> m or kg m-2]
115 
116  ! Construct source column with vanished layers removed (stored in h_nv)
117  call copy_finite_thicknesses(nz, h, cs%min_thickness, count_nonzero_layers, h_nv, mapping)
118 
119  if (count_nonzero_layers > 1) then
120  xtmp(1) = 0.0
121  do k = 1,count_nonzero_layers
122  xtmp(k+1) = xtmp(k) + h_nv(k)
123  enddo
124 
125  ! Compute densities on source column
126  pres(:) = cs%ref_pressure
127  call calculate_density(t, s, pres, densities, eqn_of_state)
128  do k = 1,count_nonzero_layers
129  densities(k) = densities(mapping(k))
130  enddo
131 
132  ! Based on source column density profile, interpolate to generate a new grid
133  call build_and_interpolate_grid(cs%interp_CS, densities, count_nonzero_layers, &
134  h_nv, xtmp, cs%target_density, cs%nk, h_new, &
135  x1, h_neglect, h_neglect_edge)
136 
137  ! Inflate vanished layers
138  call old_inflate_layers_1d(cs%min_thickness, cs%nk, h_new)
139 
140  ! Comment: The following adjustment of h_new, and re-calculation of h_new via x1 needs to be removed
141  x1(1) = 0.0 ; do k = 1,cs%nk ; x1(k+1) = x1(k) + h_new(k) ; enddo
142  do k = 1,cs%nk
143  h_new(k) = x1(k+1) - x1(k)
144  enddo
145 
146  else ! count_nonzero_layers <= 1
147  if (nz == cs%nk) then
148  h_new(:) = h(:) ! This keeps old behavior
149  else
150  h_new(:) = 0.
151  h_new(1) = h(1)
152  endif
153  endif
154 
155  ! Return interface positions
156  if (cs%integrate_downward_for_e) then
157  ! Remapping is defined integrating from zero
158  z_interface(1) = 0.
159  do k = 1,cs%nk
160  z_interface(k+1) = z_interface(k) - h_new(k)
161  enddo
162  else
163  ! The rest of the model defines grids integrating up from the bottom
164  z_interface(cs%nk+1) = -depth
165  do k = cs%nk,1,-1
166  z_interface(k) = z_interface(k+1) + h_new(k)
167  enddo
168  endif
169 
170 end subroutine build_rho_column
171 
172 !### build_rho_column_iteratively is never used or called.
173 
174 !> Iteratively build a rho coordinate column
175 !!
176 !! The algorithm operates as follows within each column:
177 !!
178 !! 1. Given T & S within each layer, the layer densities are computed.
179 !! 2. Based on these layer densities, a global density profile is reconstructed
180 !! (this profile is monotonically increasing and may be discontinuous)
181 !! 3. The new grid interfaces are determined based on the target interface
182 !! densities.
183 !! 4. T & S are remapped onto the new grid.
184 !! 5. Return to step 1 until convergence or until the maximum number of
185 !! iterations is reached, whichever comes first.
186 subroutine build_rho_column_iteratively(CS, remapCS, nz, depth, h, T, S, eqn_of_state, &
187  zInterface, h_neglect, h_neglect_edge, dev_tol)
188  type(rho_cs), intent(in) :: CS !< Regridding control structure
189  type(remapping_cs), intent(in) :: remapCS !< Remapping parameters and options
190  integer, intent(in) :: nz !< Number of levels
191  real, intent(in) :: depth !< Depth of ocean bottom [Z ~> m]
192  real, dimension(nz), intent(in) :: h !< Layer thicknesses in Z coordinates [Z ~> m]
193  real, dimension(nz), intent(in) :: T !< T for column [degC]
194  real, dimension(nz), intent(in) :: S !< S for column [ppt]
195  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
196  real, dimension(nz+1), intent(inout) :: zInterface !< Absolute positions of interfaces
197  real, optional, intent(in) :: h_neglect !< A negligibly small width for the
198  !! purpose of cell reconstructions
199  !! in the same units as h [Z ~> m]
200  real, optional, intent(in) :: h_neglect_edge !< A negligibly small width
201  !! for the purpose of edge value calculations
202  !! in the same units as h [Z ~> m]
203  real, optional, intent(in) :: dev_tol !< The tolerance for the deviation between
204  !! successive grids for determining when the
205  !! iterative solver has converged [Z ~> m]
206 
207  ! Local variables
208  real, dimension(nz+1) :: x0, x1, xTmp ! Temporary interface heights [Z ~> m]
209  real, dimension(nz) :: pres ! The pressure used in the equation of state [R L2 T-2 ~> Pa].
210  real, dimension(nz) :: densities ! Layer densities [R ~> kg m-3]
211  real, dimension(nz) :: T_tmp, S_tmp ! A temporary profile of temperature [degC] and salinity [ppt].
212  real, dimension(nz) :: Tmp ! A temporary variable holding a remapped variable.
213  real, dimension(nz) :: h0, h1, hTmp ! Temporary thicknesses [Z ~> m]
214  real :: deviation ! When iterating to determine the final grid, this is the
215  ! deviation between two successive grids [Z ~> m].
216  real :: deviation_tol ! Deviation tolerance between succesive grids in
217  ! regridding iterations [Z ~> m]
218  real :: threshold ! The minimum thickness for a layer to be considered to exist [Z ~> m]
219  integer, dimension(nz) :: mapping ! The indices of the massive layers in the initial column.
220  integer :: k, m, count_nonzero_layers
221 
222  ! Maximum number of regridding iterations
223  integer, parameter :: NB_REGRIDDING_ITERATIONS = 1
224 
225  threshold = cs%min_thickness
226  pres(:) = cs%ref_pressure
227  t_tmp(:) = t(:)
228  s_tmp(:) = s(:)
229  h0(:) = h(:)
230 
231  ! Start iterations to build grid
232  m = 1
233  deviation_tol = 1.0e-15*depth ; if (present(dev_tol)) deviation_tol = dev_tol
234 
235  do m=1,nb_regridding_iterations
236 
237  ! Construct column with vanished layers removed
238  call copy_finite_thicknesses(nz, h0, threshold, count_nonzero_layers, htmp, mapping)
239  if ( count_nonzero_layers <= 1 ) then
240  h1(:) = h0(:)
241  exit ! stop iterations here
242  endif
243 
244  xtmp(1) = 0.0
245  do k = 1,count_nonzero_layers
246  xtmp(k+1) = xtmp(k) + htmp(k)
247  enddo
248 
249  ! Compute densities within current water column
250  call calculate_density( t_tmp, s_tmp, pres, densities, eqn_of_state)
251 
252  do k = 1,count_nonzero_layers
253  densities(k) = densities(mapping(k))
254  enddo
255 
256  ! One regridding iteration
257  ! Based on global density profile, interpolate to generate a new grid
258  call build_and_interpolate_grid(cs%interp_CS, densities, count_nonzero_layers, &
259  htmp, xtmp, cs%target_density, nz, h1, x1, h_neglect, h_neglect_edge)
260 
261  call old_inflate_layers_1d( cs%min_thickness, nz, h1 )
262  x1(1) = 0.0 ; do k = 1,nz ; x1(k+1) = x1(k) + h1(k) ; enddo
263 
264  ! Remap T and S from previous grid to new grid
265  do k = 1,nz
266  h1(k) = x1(k+1) - x1(k)
267  enddo
268 
269  call remapping_core_h(remapcs, nz, h0, s, nz, h1, tmp, h_neglect, h_neglect_edge)
270  s_tmp(:) = tmp(:)
271 
272  call remapping_core_h(remapcs, nz, h0, t, nz, h1, tmp, h_neglect, h_neglect_edge)
273  t_tmp(:) = tmp(:)
274 
275  ! Compute the deviation between two successive grids
276  deviation = 0.0
277  x0(1) = 0.0
278  x1(1) = 0.0
279  do k = 2,nz
280  x0(k) = x0(k-1) + h0(k-1)
281  x1(k) = x1(k-1) + h1(k-1)
282  deviation = deviation + (x0(k)-x1(k))**2
283  enddo
284  deviation = sqrt( deviation / (nz-1) )
285 
286  if ( deviation <= deviation_tol ) exit
287 
288  ! Copy final grid onto start grid for next iteration
289  h0(:) = h1(:)
290  enddo ! end regridding iterations
291 
292  if (cs%integrate_downward_for_e) then
293  zinterface(1) = 0.
294  do k = 1,nz
295  zinterface(k+1) = zinterface(k) - h1(k)
296  ! Adjust interface position to accommodate inflating layers
297  ! without disturbing the interface above
298  enddo
299  else
300  ! The rest of the model defines grids integrating up from the bottom
301  zinterface(nz+1) = -depth
302  do k = nz,1,-1
303  zinterface(k) = zinterface(k+1) + h1(k)
304  ! Adjust interface position to accommodate inflating layers
305  ! without disturbing the interface above
306  enddo
307  endif
308 
309 end subroutine build_rho_column_iteratively
310 
311 !> Copy column thicknesses with vanished layers removed
312 subroutine copy_finite_thicknesses(nk, h_in, thresh, nout, h_out, mapping)
313  integer, intent(in) :: nk !< Number of layer for h_in, T_in, S_in
314  real, dimension(nk), intent(in) :: h_in !< Thickness of input column [H ~> m or kg m-2] or [Z ~> m]
315  real, intent(in) :: thresh !< Thickness threshold defining vanished
316  !! layers [H ~> m or kg m-2] or [Z ~> m]
317  integer, intent(out) :: nout !< Number of non-vanished layers
318  real, dimension(nk), intent(out) :: h_out !< Thickness of output column [H ~> m or kg m-2] or [Z ~> m]
319  integer, dimension(nk), intent(out) :: mapping !< Index of k-out corresponding to k-in
320  ! Local variables
321  integer :: k, k_thickest
322  real :: thickness_in_vanished ! Summed thicknesses in discarded layers [H ~> m or kg m-2] or [Z ~> m]
323  real :: thickest_h_out ! Thickness of the thickest layer [H ~> m or kg m-2] or [Z ~> m]
324 
325  ! Build up new grid
326  nout = 0
327  thickness_in_vanished = 0.0
328  thickest_h_out = h_in(1)
329  k_thickest = 1
330  do k = 1, nk
331  mapping(k) = nout ! Note k>=nout always
332  h_out(k) = 0. ! Make sure h_out is set everywhere
333  if (h_in(k) > thresh) then
334  ! For non-vanished layers
335  nout = nout + 1
336  mapping(nout) = k
337  h_out(nout) = h_in(k)
338  if (h_out(nout) > thickest_h_out) then
339  thickest_h_out = h_out(nout)
340  k_thickest = nout
341  endif
342  else
343  ! Add up mass in vanished layers
344  thickness_in_vanished = thickness_in_vanished + h_in(k)
345  endif
346  enddo
347 
348  ! No finite layers
349  if (nout <= 1) return
350 
351  ! Adjust for any lost volume in vanished layers
352  h_out(k_thickest) = h_out(k_thickest) + thickness_in_vanished
353 
354 end subroutine copy_finite_thicknesses
355 
356 !------------------------------------------------------------------------------
357 !> Inflate vanished layers to finite (nonzero) width
358 subroutine old_inflate_layers_1d( min_thickness, nk, h )
360  ! Argument
361  real, intent(in) :: min_thickness !< Minimum allowed thickness [H ~> m or kg m-2]
362  integer, intent(in) :: nk !< Number of layers in the grid
363  real, dimension(:), intent(inout) :: h !< Layer thicknesses [H ~> m or kg m-2]
364 
365  ! Local variable
366  integer :: k
367  integer :: k_found
368  integer :: count_nonzero_layers
369  real :: delta
370  real :: correction
371  real :: maxThickness
372 
373  ! Count number of nonzero layers
374  count_nonzero_layers = 0
375  do k = 1,nk
376  if ( h(k) > min_thickness ) then
377  count_nonzero_layers = count_nonzero_layers + 1
378  endif
379  enddo
380 
381  ! If all layer thicknesses are greater than the threshold, exit routine
382  if ( count_nonzero_layers == nk ) return
383 
384  ! If all thicknesses are zero, inflate them all and exit
385  if ( count_nonzero_layers == 0 ) then
386  do k = 1,nk
387  h(k) = min_thickness
388  enddo
389  return
390  endif
391 
392  ! Inflate zero layers
393  correction = 0.0
394  do k = 1,nk
395  if ( h(k) <= min_thickness ) then
396  delta = min_thickness - h(k)
397  correction = correction + delta
398  h(k) = h(k) + delta
399  endif
400  enddo
401 
402  ! Modify thicknesses of nonzero layers to ensure volume conservation
403  maxthickness = h(1)
404  k_found = 1
405  do k = 1,nk
406  if ( h(k) > maxthickness ) then
407  maxthickness = h(k)
408  k_found = k
409  endif
410  enddo
411 
412  h(k_found) = h(k_found) - correction
413 
414 end subroutine old_inflate_layers_1d
415 
416 end module coord_rho
Control structure containing required parameters for the rho coordinate.
Definition: coord_rho.F90:14
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68
Provides column-wise vertical remapping functions.
Container for remapping parameters.
Routines for error handling and I/O management.
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
Control structure for regrid_interp module.
Regrid columns for the continuous isopycnal (rho) coordinate.
Definition: coord_rho.F90:2
Vertical interpolation for regridding.
A control structure for the equation of state.
Definition: MOM_EOS.F90:108