check_inputs Subroutine

public pure subroutine check_inputs(n, m, np, q, isodr, anajac, beta, ifixb, ldifx, ldscld, ldstpd, ldwe, ld2we, ldwd, ld2wd, lrwork, lrwkmin, liwork, liwkmin, sclb, scld, stpb, stpd, info, lower, upper)

Uses

  • proc~~check_inputs~~UsesGraph proc~check_inputs check_inputs module~odrpack_kinds odrpack_kinds proc~check_inputs->module~odrpack_kinds iso_fortran_env iso_fortran_env module~odrpack_kinds->iso_fortran_env

Check input parameters, indicating errors found using nonzero values of argument info.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

Number of observations.

integer, intent(in) :: m

Number of columns of data in the explanatory variable.

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

logical, intent(in) :: isodr

Variable designating whether the solution is by ODR (.true.) or by OLS (.false.).

logical, intent(in) :: anajac

Variable designating whether the Jacobians are computed by finite differences (.false.) or not (.true.).

real(kind=wp), intent(in) :: beta(np)

Function parameters.

integer, intent(in) :: ifixb(np)

Values designating whether the elements of beta are fixed at their input values or not.

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

integer, intent(in) :: ldscld

Leading dimension of array scld.

integer, intent(in) :: ldstpd

Leading dimension of array stpd.

integer, intent(in) :: ldwe

Leading dimension of array we.

integer, intent(in) :: ld2we

Second dimension of array we.

integer, intent(in) :: ldwd

Leading dimension of array wd.

integer, intent(in) :: ld2wd

Second dimension of array wd.

integer, intent(in) :: lrwork

Length of vector rwork.

integer, intent(in) :: lrwkmin

Minimum acceptable length of array rwork.

integer, intent(in) :: liwork

Length of vector iwork.

integer, intent(in) :: liwkmin

Minimum acceptable length of array iwork.

real(kind=wp), intent(in) :: sclb(np)

Scaling values for beta.

real(kind=wp), intent(in) :: scld(ldscld,m)

Scaling value for delta.

real(kind=wp), intent(in) :: stpb(np)

Step for the finite difference derivative wrt beta.

real(kind=wp), intent(in) :: stpd(ldstpd,m)

Step for the finite difference derivative wrt delta.

integer, intent(out) :: info

Variable designating why the computations were stopped.

real(kind=wp), intent(in) :: lower(np)

Lower bound on beta.

real(kind=wp), intent(in) :: upper(np)

Upper bound on beta.


Called by

proc~~check_inputs~~CalledByGraph proc~check_inputs check_inputs proc~odr odr proc~odr->proc~check_inputs proc~odr_long_c odr_long_c proc~odr_long_c->proc~odr proc~odr_medium_c odr_medium_c proc~odr_medium_c->proc~odr proc~odr_short_c odr_short_c proc~odr_short_c->proc~odr program~example1 example1 program~example1->proc~odr program~example2 example2 program~example2->proc~odr program~example3 example3 program~example3->proc~odr program~example4 example4 program~example4->proc~odr program~example5 example5 program~example5->proc~odr

Variables

Type Visibility Attributes Name Initial
integer, public :: last
integer, public :: npp

Source Code

   pure subroutine check_inputs &
      (n, m, np, q, &
       isodr, anajac, &
       beta, ifixb, &
       ldifx, ldscld, ldstpd, ldwe, ld2we, ldwd, ld2wd, &
       lrwork, lrwkmin, liwork, liwkmin, &
       sclb, scld, stpb, stpd, &
       info, &
       lower, upper)
   !! Check input parameters, indicating errors found using nonzero values of argument `info`.

      use odrpack_kinds, only: zero

      integer, intent(in) :: n
         !! Number of observations.
      integer, intent(in) :: m
         !! Number of columns of data in the explanatory variable.
      integer, intent(in) :: np
         !! Number of function parameters.
      integer, intent(in) :: q
         !! Number of responses per observation.
      logical, intent(in) :: isodr
         !! Variable designating whether the solution is by ODR (`.true.`) or by OLS (`.false.`).
      logical, intent(in) :: anajac
         !! Variable designating whether the Jacobians are computed by finite differences
         !! (`.false.`) or not (`.true.`).
      real(wp), intent(in) :: beta(np)
         !! Function parameters.
      integer, intent(in) :: ifixb(np)
         !! Values designating whether the elements of `beta` are fixed at their input values or not.
      integer, intent(in) :: ldifx
         !! Leading dimension of array `ifixx`.
      integer, intent(in) :: ldscld
         !! Leading dimension of array `scld`.
      integer, intent(in) :: ldstpd
         !! Leading dimension of array `stpd`.
      integer, intent(in) :: ldwe
      !! Leading dimension of array `we`.
      integer, intent(in) :: ld2we
         !! Second dimension of array `we`.
      integer, intent(in) :: ldwd
         !! Leading dimension of array `wd`.
      integer, intent(in) :: ld2wd
         !! Second dimension of array `wd`.
      integer, intent(in) :: lrwork
         !! Length of vector `rwork`.
      integer, intent(in) :: lrwkmin
         !! Minimum acceptable length of array `rwork`.
      integer, intent(in) :: liwork
         !! Length of vector `iwork`.
      integer, intent(in) :: liwkmin
         !! Minimum acceptable length of array `iwork`.
      real(wp), intent(in) :: sclb(np)
         !! Scaling values for `beta`.
      real(wp), intent(in) :: scld(ldscld, m)
         !! Scaling value for `delta`.
      real(wp), intent(in) :: stpb(np)
         !! Step for the finite difference derivative wrt `beta`.
      real(wp), intent(in) :: stpd(ldstpd, m)
         !! Step for the finite difference derivative wrt `delta`.
      integer, intent(out) :: info
         !! Variable designating why the computations were stopped.
      real(wp), intent(in) :: lower(np)
         !! Lower bound on `beta`.
      real(wp), intent(in) :: upper(np)
         !! Upper bound on `beta`.

      ! Local scalars
      integer :: last, npp

      ! Variable Definitions (alphabetically)
      !  LAST:    The last row of the array to be accessed.
      !  NPP:     The number of function parameters being estimated.

      ! Find actual number of parameters being estimated
      if ((np <= 0) .or. (ifixb(1) < 0)) then
         npp = np
      else
         npp = count(ifixb /= 0)
      end if

      ! Check problem specification parameters
      if ((n <= 0) .or. (m <= 0) .or. (npp <= 0 .or. npp > n) .or. (q <= 0)) then
         info = 10000
         if (n <= 0) then
            info = info + 1000
         end if
         if (m <= 0) then
            info = info + 100
         end if
         if (npp <= 0 .or. npp > n) then
            info = info + 10
         end if
         if (q <= 0) then
            info = info + 1
         end if
         return
      end if

      ! Check dimension specification parameters
      if (((ldwe /= 1) .and. (ldwe < n)) .or. &
          ((ld2we /= 1) .and. (ld2we < q)) .or. &
          (isodr .and. ((ldwd /= 1) .and. (ldwd < n))) .or. &
          (isodr .and. ((ld2wd /= 1) .and. (ld2wd < m))) .or. &
          (isodr .and. ((ldifx /= 1) .and. (ldifx < n))) .or. &
          (isodr .and. ((ldstpd /= 1) .and. (ldstpd < n))) .or. &
          (isodr .and. ((ldscld /= 1) .and. (ldscld < n))) .or. &
          (lrwork < lrwkmin) .or. &
          (liwork < liwkmin)) then

         info = 20000

         if ((ldwe /= 1 .and. ldwe < n) .or. (ld2we /= 1 .and. ld2we < q)) then
            info = info + 100
         end if

         if (isodr .and. &
             ((ldwd /= 1 .and. ldwd < n) .or. (ld2wd /= 1 .and. ld2wd < m))) then
            info = info + 200
         end if

         if (isodr .and. (ldifx /= 1 .and. ldifx < n)) then
            info = info + 10
         end if

         if (isodr .and. (ldstpd /= 1 .and. ldstpd < n)) then
            info = info + 20
         end if

         if (isodr .and. &
             (ldscld /= 1 .and. ldscld < n)) then
            info = info + 40
         end if

         if (lrwork < lrwkmin) then
            info = info + 1
         end if

         if (liwork < liwkmin) then
            info = info + 2
         end if

      end if

      ! Check DELTA scaling
      if (isodr .and. scld(1, 1) > zero) then
         if (ldscld >= n) then
            last = n
         else
            last = 1
         end if
         if (any(scld(1:last, :) <= zero)) then
            info = 30200
         end if
      end if

      ! Check BETA scaling
      if (sclb(1) > zero) then
         if (any(sclb <= zero)) then
            if (info == 0) then
               info = 30100
            else
               info = info + 100
            end if
         end if
      end if

      ! Check DELTA finite difference step sizes
      if ((.not. anajac) .and. isodr .and. stpd(1, 1) > zero) then
         if (ldstpd >= n) then
            last = n
         else
            last = 1
         end if
         if (any(stpd(1:last, :) <= zero)) then
            if (info == 0) then
               info = 32000
            else
               info = info + 2000
            end if
         end if
      end if

      ! Check BETA finite difference step sizes
      if ((.not. anajac) .and. stpb(1) > zero) then
         if (any(stpb <= zero)) then
            if (info == 0) then
               info = 31000
            else
               info = info + 1000
            end if
         end if
      end if

      !  Check bounds
      if (any(upper < lower)) then
         if (info == 0) then
            info = 91000
         end if
      end if

      if (any( &
          ((upper < beta) .or. (lower > beta)) &
          .and. .not. (upper < lower))) then
         if (info >= 90000) then
            info = info + 100
         else
            info = 90100
         end if
      end if

   end subroutine check_inputs