Skip to content
Snippets Groups Projects
get_1_outerm.f90 7.85 KiB
Newer Older
module get_1_outerm_m
  subroutine get_1_outerm(out_cont, cont_list, cont_list_proj, cont_list_ssh, &
       n_cont, cyclone, extr_coord_proj, innermost_level_2, outside_points, &

    ! This procedure gets one outermost good contour with sufficient
    ! area, if it can find one. If the procedure cannot find an
    ! outermost contour, it returns a null ssh_contour. The procedure
Lionel GUEZ's avatar
Lionel GUEZ committed
    ! also returns some inner good contours. The maximum number of
    ! inner good contours returned is given by the size of cont_list -
    ! 1. If we find more inner good contours, then we return only the
    ! first ones found, because they provide a better sample of the
    ! ssh range. If the procedure cannot find an outermost contour, it
    ! does not return any inner contour. An outermost contour is not
    ! returned if there is no good contour at innermost_level_2 or if
    ! the outermost contour does not have sufficient area. If an
    ! outermost contour is found then its level is farther from the
    ! extremum than innermost_level_2. The outermost contour is
GUEZ Lionel's avatar
GUEZ Lionel committed
    ! duplicated in the array cont_list. In cont_list_proj, the
    ! contours, including the outermost contour, are in monotonic
    ! order of ssh: ascending if the extremum is a minimum, descending
    ! if the extremum is a maximum. On return, the array
    ! abs(ediff1d(cont_list_ssh(:n_cont))) is in descending
Lionel GUEZ's avatar
Lionel GUEZ committed
    ! order. (See notes for a proof.) On return, n_cont == 1 means we
    ! have only found one good contour at innermost_level_2 and it has
GUEZ Lionel's avatar
GUEZ Lionel committed
    ! sufficient area. On return, n_cont == 0 .eqv. out_cont ==
    ! null_polyline().
Lionel GUEZ's avatar
Lionel GUEZ committed

    ! Method: dichotomy on level of ssh.
GUEZ Lionel's avatar
GUEZ Lionel committed
    ! The length of the interval of longitudes of the domain,
    ! longitude of ssh(size(ssh, 1), 1) - longitude of ssh(1, 1),
    ! should be < 180 degrees, so that the geometrical processing done
    ! here is non-ambiguous. The longitudes of outside points and the
    ! target extremum should be in the interval [longitude of ssh(1,
    ! 1), longitude of ssh(size(ssh, 1), 1)] as there will be no
    ! automatic shifting by a mulitple of 360 degrees.
    use contour_531, only: polyline, convert_to_reg_coord
    use jumble, only: assert
Lionel GUEZ's avatar
Lionel GUEZ committed
    use ccw_orient_m, only: ccw_orient
    use config_m, only: min_area
Lionel GUEZ's avatar
Lionel GUEZ committed
    use derived_types, only: ssh_contour, null_ssh_contour
    use good_contour_m, only: good_contour
    use input_ssh_m, only: corner_whole, step
    use spher_polyline_area_m, only: spher_polyline_area
    type(ssh_contour), intent(out):: out_cont
Lionel GUEZ's avatar
Lionel GUEZ committed
    ! points given by longitude, latitude, in rad
GUEZ Lionel's avatar
GUEZ Lionel committed
    type(polyline), intent(out):: cont_list(:) ! (n_max_cont)
Lionel GUEZ's avatar
Lionel GUEZ committed
    ! Contour list, with points given by longitude and latitude, in
    ! radians. Defined only for subscripts 1:n_cont. If n_cont >= 1,
GUEZ Lionel's avatar
GUEZ Lionel committed
    ! the outermost contour is element with subscript n_cont. For i <=
    ! n_cont, the order of points in cont_list(i) (clockwise or
Lionel GUEZ's avatar
Lionel GUEZ committed
    ! counter-clockwise) is not specified.
Lionel GUEZ's avatar
Lionel GUEZ committed

    type(polyline), intent(out):: cont_list_proj(:) ! (n_max_cont)
    ! Contour list, with points given by projection
    ! coordinates. Defined only for subscripts 1:n_cont. If n_cont >=
    ! 1, the outermost contour is element with subscript n_cont. The
GUEZ Lionel's avatar
GUEZ Lionel committed
    ! contours are in monotonic order of ssh. The contour with a given
    ! subscript in cont_list_proj corresponds to the contour with the
GUEZ Lionel's avatar
GUEZ Lionel committed
    ! same subscript in cont_list. More precisely, the point with
    ! subscript j in contour i of cont_list_proj corresponds to the
    ! point with subscript j in contour i of cont_list (this is
    ! important for mean_speed). For i <= n_cont, the order of points
    ! in cont_list_proj(i) (clockwise or counter-clockwise) is not
    ! specified.

    real, intent(out):: cont_list_ssh(:) ! (n_max_cont)
Lionel GUEZ's avatar
Lionel GUEZ committed
    integer, intent(out):: n_cont
Lionel GUEZ's avatar
Lionel GUEZ committed
    ! number of good contours found and stored, 0 <= n_cont <= n_max_cont - 1

    real, intent(in):: extr_coord_proj(:) ! (2)
    ! coordinates of extremum in projection space
    real, intent(in):: innermost_level_2
    ! ssh level of the innermost contour to consider around the target
    ! extremum, in m. Assume: innermost_level_2 < extremum for a
    ! maximum, > extremum for a minimum.
    real, intent(in):: outside_points(:, :) ! (2, :)
    ! coordinates in projection space of all the significant extrema,
    ! except the target extremum
    real, intent(in):: ssh(:, :) ! in m

    real, intent(in):: corner(:) ! (2)
    ! Coordinates in projection space corresponding to ssh(1,1). A
Lionel GUEZ's avatar
Lionel GUEZ committed
    ! priori, this is not the corner of the grid of the NetCDF file.

    ! Local:
    real level_try, level_good, level_bad ! in m
    real, parameter:: precision = 5e-4 ! in m
    logical mask(size(ssh, 1), size(ssh, 2))
Lionel GUEZ's avatar
Lionel GUEZ committed
    integer n_max_cont ! >= 3

    !-----------------------------------------------------------------

GUEZ Lionel's avatar
GUEZ Lionel committed
    n_max_cont = size(cont_list_proj)
Lionel GUEZ's avatar
Lionel GUEZ committed
    n_cont = 0
    cont_list_proj(1) = good_contour(corner, ssh, innermost_level_2, &
         extr_coord_proj, outside_points)
    test_n_points: if (.not. cont_list_proj(1)%closed) then
GUEZ Lionel's avatar
GUEZ Lionel committed
       ! cont_list_proj(1)%n_points == 0
       out_cont = null_ssh_contour()
Lionel GUEZ's avatar
Lionel GUEZ committed
    else test_n_points
Lionel GUEZ's avatar
Lionel GUEZ committed
       n_cont = 1
       cont_list_ssh(1) = innermost_level_2
       level_good = innermost_level_2
       mask = ssh /= huge(0.)
       level_try = merge(maxval(ssh, mask), minval(ssh, mask), cyclone)
       call assert(merge(level_try > level_good, level_try < level_good, &
            cyclone), "get_1_outerm level_try")
       cont_list_proj(2) = good_contour(corner, ssh, level_try, &
            extr_coord_proj, outside_points)
       if (cont_list_proj(2)%closed) then
GUEZ Lionel's avatar
GUEZ Lionel committed
          ! cont_list_proj(2)%n_points /= 0
          ! This can only happen if there are degenerate extrema
          ! everywhere around the path of contour. So it should
          ! probably never happen.
Lionel GUEZ's avatar
Lionel GUEZ committed
          n_cont = 2
          cont_list_ssh(2) = level_try
Lionel GUEZ's avatar
Lionel GUEZ committed
          ! Assertion: {no good contour at maxval or minval}
          level_bad = level_try

          do while (abs(level_bad - level_good) > precision)
             level_try = (level_good + level_bad) / 2.
             cont_list_proj(n_cont + 1) = good_contour(corner, ssh, level_try, &
                  extr_coord_proj, outside_points)
             if (cont_list_proj(n_cont + 1)%closed) then
GUEZ Lionel's avatar
GUEZ Lionel committed
                ! cont_list_proj(n_cont + 1)%n_points /= 0
                cont_list_ssh(n_cont + 1) = level_try
Lionel GUEZ's avatar
Lionel GUEZ committed
                if (n_cont == n_max_cont - 1) then
Lionel GUEZ's avatar
Lionel GUEZ committed
                   ! Replace the previous good contour found
                   cont_list_proj(n_cont) = cont_list_proj(n_max_cont)
Lionel GUEZ's avatar
Lionel GUEZ committed
                else
                   n_cont = n_cont + 1
                end if
Lionel GUEZ's avatar
Lionel GUEZ committed

                level_good = level_try
             else
                level_bad = level_try
             end if

             ! Assertion: cont_list_proj(n_cont) ==
             ! good_contour(. . ., level_good) and
             ! cont_list_proj(n_cont)%closed and no good contour at
             ! level_bad
Lionel GUEZ's avatar
Lionel GUEZ committed
          ! Assertion: n_cont <= n_max_cont - 1
GUEZ Lionel's avatar
GUEZ Lionel committed
       out_cont%polyline &
            = convert_to_reg_coord(cont_list_proj(n_cont), corner_whole, step)
       out_cont%area = spher_polyline_area(out_cont%polyline)
       if (abs(out_cont%area) < min_area) then
          out_cont = null_ssh_contour()
          n_cont = 0
       else
          ! n_cont >= 1
GUEZ Lionel's avatar
GUEZ Lionel committed
          forall (i = 1:n_cont - 1) cont_list(i) &
               = convert_to_reg_coord(cont_list_proj(i), corner_whole, step)
GUEZ Lionel's avatar
GUEZ Lionel committed
          cont_list(n_cont) = out_cont%polyline
          out_cont%ssh = cont_list_ssh(n_cont)
GUEZ Lionel's avatar
GUEZ Lionel committed

          call ccw_orient(out_cont)
GUEZ Lionel's avatar
GUEZ Lionel committed
          ! (Note we have to orient out_cont after copying it to
          ! cont_list, so that cont_list corresponds to cont_list_proj
          ! point by point.)
Lionel GUEZ's avatar
Lionel GUEZ committed
    end if test_n_points
  end subroutine get_1_outerm
end module get_1_outerm_m