dodchk Subroutine

public pure subroutine dodchk(n, m, np, nq, isodr, anajac, implct, beta, ifixb, ldx, ldifx, ldscld, ldstpd, ldwe, ld2we, ldwd, ld2wd, ldy, lwork, lwkmn, liwork, liwkmn, sclb, scld, stpb, stpd, info, lower, upper)

Uses

  • proc~~dodchk~~UsesGraph proc~dodchk dodchk module~odrpack_kinds odrpack_kinds proc~dodchk->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

The number of observations.

integer, intent(in) :: m

The number of columns of data in the explanatory variable.

integer, intent(in) :: np

The number of function parameters.

integer, intent(in) :: nq

The number of responses per observation.

logical, intent(in) :: isodr

The variable designating whether the solution is by ODR (isodr = .true.) or by OLS (isodr = .false.).

logical, intent(in) :: anajac

The variable designating whether the Jacobians are computed by finite differences (anajac = .false.) or not (anajac = .true.).

logical, intent(in) :: implct

The variable designating whether the solution is by implicit ODR (implct = .true.) or explicit ODR (implct = .false.).

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

The function parameters.

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

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

integer, intent(in) :: ldx

The leading dimension of array x.

integer, intent(in) :: ldifx

The leading dimension of array ifixx.

integer, intent(in) :: ldscld

The leading dimension of array scld.

integer, intent(in) :: ldstpd

The leading dimension of array stpd.

integer, intent(in) :: ldwe

The leading dimension of array we.

integer, intent(in) :: ld2we

The second dimension of array we.

integer, intent(in) :: ldwd

The leading dimension of array wd.

integer, intent(in) :: ld2wd

The second dimension of array wd.

integer, intent(in) :: ldy

The leading dimension of array y.

integer, intent(in) :: lwork

The length of vector work.

integer, intent(in) :: lwkmn

The minimum acceptable length of array work.

integer, intent(in) :: liwork

The length of vector iwork.

integer, intent(in) :: liwkmn

The minimum acceptable length of array iwork.

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

The scaling values for beta.

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

The scaling value for delta.

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

The step for the finite difference derivative wrt beta.

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

The step for the finite difference derivative wrt delta.

integer, intent(out) :: info

The variable designating why the computations were stopped.

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

The lower bound on beta.

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

The upper bound on beta.


Called by

proc~~dodchk~~CalledByGraph proc~dodchk dodchk proc~doddrv doddrv proc~doddrv->proc~dodchk proc~dodcnt dodcnt proc~dodcnt->proc~doddrv proc~odr odr proc~odr->proc~dodcnt 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 dodchk &
      (n, m, np, nq, &
       isodr, anajac, implct, &
       beta, ifixb, &
       ldx, ldifx, ldscld, ldstpd, ldwe, ld2we, ldwd, ld2wd, &
       ldy, &
       lwork, lwkmn, liwork, liwkmn, &
       sclb, scld, stpb, stpd, &
       info, &
       lower, upper)
   !! Check input parameters, indicating errors found using nonzero values of argument `info`.
      ! Routines Called  (None)
      ! Date Written   860529   (YYMMDD)
      ! Revision Date  920619   (YYMMDD)

      use odrpack_kinds, only: zero

      integer, intent(in) :: n
         !! The number of observations.
      integer, intent(in) :: m
         !! The number of columns of data in the explanatory variable.
      integer, intent(in) :: np
         !! The number of function parameters.
      integer, intent(in) :: nq
         !! The number of responses per observation.
      logical, intent(in) :: isodr
         !! The variable designating whether the solution is by ODR (`isodr = .true.`) or
         !! by OLS (`isodr = .false.`).
      logical, intent(in) :: anajac
         !! The variable designating whether the Jacobians are computed by finite differences
         !! (`anajac = .false.`) or not (`anajac = .true.`).
      logical, intent(in) :: implct
         !! The variable designating whether the solution is by implicit ODR (`implct = .true.`)
         !! or explicit ODR (`implct = .false.`).
      real(wp), intent(in) :: beta(np)
         !! The function parameters.
      integer, intent(in) :: ifixb(np)
         !! The values designating whether the elements of `beta` are fixed at their input values or not.
      integer, intent(in) :: ldx
         !! The leading dimension of array `x`.
      integer, intent(in) :: ldifx
         !! The leading dimension of array `ifixx`.
      integer, intent(in) :: ldscld
         !! The leading dimension of array `scld`.
      integer, intent(in) :: ldstpd
         !! The leading dimension of array `stpd`.
      integer, intent(in) :: ldwe
      !! The leading dimension of array `we`.
      integer, intent(in) :: ld2we
         !! The second dimension of array `we`.
      integer, intent(in) :: ldwd
         !! The leading dimension of array `wd`.
      integer, intent(in) :: ld2wd
         !! The second dimension of array `wd`.
      integer, intent(in) :: ldy
         !! The leading dimension of array `y`.
      integer, intent(in) :: lwork
         !! The length of vector `work`.
      integer, intent(in) :: lwkmn
         !! The minimum acceptable length of array `work`.
      integer, intent(in) :: liwork
         !! The length of vector `iwork`.
      integer, intent(in) :: liwkmn
         !! The minimum acceptable length of array `iwork`.
      real(wp), intent(in) :: sclb(np)
         !! The scaling values for `beta`.
      real(wp), intent(in) :: scld(ldscld, m)
         !! The scaling value for `delta`.
      real(wp), intent(in) :: stpb(np)
         !! The step for the finite difference derivative wrt `beta`.
      real(wp), intent(in) :: stpd(ldstpd, m)
         !! The step for the finite difference derivative wrt `delta`.
      integer, intent(out) :: info
         !! The variable designating why the computations were stopped.
      real(wp), intent(in) :: lower(np)
         !! The lower bound on `beta`.
      real(wp), intent(in) :: upper(np)
         !! The upper bound on `beta`.

      ! Local scalars
      integer :: last, npp

      ! Variable Definitions (alphabetically)
      !  ANAJAC:  The variable designating whether the Jacobians are computed by finite
      !           differences (ANAJAC=FALSE) or not (ANAJAC=TRUE).
      !  I:       An indexing variable.
      !  IFIXB:   The values designating whether the elements of BETA are fixed at their input
      !           values or not.
      !  IMPLCT:  The variable designating whether the solution is by implicit ODR (IMPLCT=TRUE)
      !           or explicit ODR (IMPLCT=FALSE).
      !  INFO:    The variable designating why the computations were stopped.
      !  ISODR:   The variable designating whether the solution is by ODR (ISODR=TRUE) or
      !           by OLS (ISODR=FALSE).
      !  J:       An indexing variable.
      !  K:       An indexing variable.
      !  LAST:    The last row of the array to be accessed.
      !  LDIFX:   The leading dimension of array IFIXX.
      !  LDSCLD:  The leading dimension of array SCLD.
      !  LDSTPD:  The leading dimension of array STPD.
      !  LDWD:    The leading dimension of array WD.
      !  LDWE:    The leading dimension of array WE.
      !  LDX:     The leading dimension of array X.
      !  LDY:     The leading dimension of array X.
      !  LD2WD:   The second dimension of array WD.
      !  LD2WE:   The second dimension of array WE.
      !  LIWKMN:  The minimum acceptable length of array IWORK.
      !  LIWORK:  The length of vector IWORK.
      !  LWKMN:   The minimum acceptable length of array WORK.
      !  LWORK:   The length of vector WORK.
      !  M:       The number of columns of data in the explanatory variable.
      !  N:       The number of observations.
      !  NP:      The number of function parameters.
      !  NPP:     The number of function parameters being estimated.
      !  NQ:      The number of responses per observations.
      !  SCLB:    The scaling values for BETA.
      !  SCLD:    The scaling value for DELTA.
      !  STPB:    The step for the finite difference derivitive wrt BETA.
      !  STPD:    The step for the finite difference derivitive wrt DELTA.

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

      ! Check problem specification parameters
      if ((n <= 0) .or. (m <= 0) .or. (npp <= 0 .or. npp > n) .or. (nq <= 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 (nq <= 0) then
            info = info + 1
         end if
         return
      end if

      ! Check dimension specification parameters
      if ((.not. implct .and. (ldy < n)) .or. &
          (ldx < n) .or. &
          ((ldwe /= 1) .and. (ldwe < n)) .or. &
          ((ld2we /= 1) .and. (ld2we < nq)) .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. &
          (lwork < lwkmn) .or. &
          (liwork < liwkmn)) then

         info = 20000

         if (.not. implct .and. ldy < n) then
            info = info + 1000
         end if

         if (ldx < n) then
            info = info + 2000
         end if

         if ((ldwe /= 1 .and. ldwe < n) .or. (ld2we /= 1 .and. ld2we < nq)) 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 (lwork < lwkmn) then
            info = info + 1
         end if

         if (liwork < liwkmn) 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, 1:m) <= zero)) then
            info = 30200
         end if
      end if

      ! Check BETA scaling
      if (sclb(1) > zero) then
         if (any(sclb(1:np) <= 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 (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, 1:m) <= 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 (anajac .and. stpb(1) > zero) then
         if (any(stpb(1:np) <= zero)) then
            if (info == 0) then
               info = 31000
            else
               info = info + 1000
            end if
         end if
      end if

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

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

   end subroutine dodchk