dodpcr Subroutine

public impure subroutine dodpcr(ipr, lunrpt, head, prtpen, fstitr, didvcv, iflag, n, m, np, nq, npp, nnzw, msgb, msgd, beta, y, ldy, x, ldx, delta, we, ldwe, ld2we, wd, ldwd, ld2wd, ifixb, ifixx, ldifx, lower, upper, ssf, tt, ldtt, stpb, stpd, ldstpd, job, neta, taufac, sstol, partol, maxit, wss, rvar, idf, sdbeta, niter, nfev, njev, actred, prered, tau, pnorm, alpha, f, rcond, irank, info, istop)

Uses

  • proc~~dodpcr~~UsesGraph proc~dodpcr dodpcr module~odrpack_core odrpack_core proc~dodpcr->module~odrpack_core module~odrpack_kinds odrpack_kinds module~odrpack_core->module~odrpack_kinds iso_fortran_env iso_fortran_env module~odrpack_kinds->iso_fortran_env

Generate computation reports.

Arguments

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

The variable indicating what is to be printed.

integer, intent(in) :: lunrpt

The logical unit number used for computation reports.

logical, intent(inout) :: head

The variable designating whether the heading is to be printed (head = .true.) or not (head = .false.).

logical, intent(in) :: prtpen

The variable designating whether the penalty parameter is to be printed in the iteration report (prtpen = .true.) or not (prtpen = .false.).

logical, intent(in) :: fstitr

The variable designating whether this is the first iteration (fstitr = .true.) or not (fstitr = .false.).

logical, intent(in) :: didvcv

The variable designating whether the covariance matrix was computed (didvcv = .true.) or not (didvcv = .false.).

integer, intent(in) :: iflag

The variable designating what is to be printed.

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.

integer, intent(in) :: npp

The number of function parameters being estimated.

integer, intent(in) :: nnzw

The number of nonzero weighted observations.

integer, intent(in) :: msgb(nq*np+1)

The error checking results for the Jacobian with respect to beta.

integer, intent(in) :: msgd(nq*m+1)

The error checking results for the Jacobian with respect to delta.

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

The function parameters.

real(kind=wp), intent(in) :: y(ldy,nq)

The response variable. Unused when the model is implicit.

integer, intent(in) :: ldy

The leading dimension of array y.

real(kind=wp), intent(in) :: x(ldx,m)

The explanatory variable.

integer, intent(in) :: ldx

The leading dimension of array x.

real(kind=wp), intent(in) :: delta(n,m)

The estimated errors in the explanatory variables.

real(kind=wp), intent(in) :: we(ldwe,ld2we,nq)

The epsilon weights.

integer, intent(in) :: ldwe

The leading dimension of array we.

integer, intent(in) :: ld2we

The second dimension of array we.

real(kind=wp), intent(in) :: wd(ldwd,ld2wd,m)

The delta weights.

integer, intent(in) :: ldwd

The leading dimension of array wd.

integer, intent(in) :: ld2wd

The second dimension of array wd.

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

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

integer, intent(in) :: ifixx(ldifx,m)

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

integer, intent(in) :: ldifx

The leading dimension of array ifixx.

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

The lower bounds for beta.

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

The upper bounds for beta.

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

The scaling values for beta.

real(kind=wp), intent(in) :: tt(ldtt,m)

The scaling values for delta.

integer, intent(in) :: ldtt

The leading dimension of array tt.

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

The relative step used for computing finite difference derivatives with respect to beta.

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

The relative step used for computing finite difference derivatives with respect to delta.

integer, intent(in) :: ldstpd

The leading dimension of array stpd.

integer, intent(in) :: job

The variable controlling problem initialization and computational method.

integer, intent(in) :: neta

The number of accurate digits in the function results.

real(kind=wp), intent(in) :: taufac

The factor used to compute the initial trust region diameter.

real(kind=wp), intent(in) :: sstol

The sum-of-squares convergence stopping tolerance.

real(kind=wp), intent(in) :: partol

The parameter convergence stopping tolerance.

integer, intent(in) :: maxit

The maximum number of iterations allowed.

real(kind=wp), intent(in) :: wss(3)

The sum-of-squares of the weighted epsilons and deltas, the sum-of-squares of the weighted deltas, and the sum-of-squares of the weighted epsilons.

real(kind=wp), intent(in) :: rvar

The residual variance.

integer, intent(in) :: idf

The degrees of freedom of the fit, equal to the number of observations with nonzero weighted derivatives minus the number of parameters being estimated.

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

The standard errors of the estimated parameters.

integer, intent(in) :: niter

The number of iterations.

integer, intent(in) :: nfev

The number of function evaluations.

integer, intent(in) :: njev

The number of Jacobian evaluations.

real(kind=wp), intent(in) :: actred

The actual relative reduction in the sum-of-squares.

real(kind=wp), intent(in) :: prered

The predicted relative reduction in the sum-of-squares.

real(kind=wp), intent(in) :: tau

The trust region diameter.

real(kind=wp), intent(in) :: pnorm

The norm of the scaled estimated parameters.

real(kind=wp), intent(in) :: alpha

The Levenberg-Marquardt parameter.

real(kind=wp), intent(in) :: f(n,nq)

The estimated values of epsilon.

real(kind=wp), intent(in) :: rcond

The approximate reciprocal condition of tfjacb.

integer, intent(in) :: irank

The rank deficiency of the Jacobian with respect to beta.

integer, intent(in) :: info

The variable designating why the computations were stopped.

integer, intent(in) :: istop

The variable designating whether there are problems computing the function at the current beta and delta.


Calls

proc~~dodpcr~~CallsGraph proc~dodpcr dodpcr proc~dflags dflags proc~dodpcr->proc~dflags proc~dodpc1 dodpc1 proc~dodpcr->proc~dodpc1 proc~dodpc2 dodpc2 proc~dodpcr->proc~dodpc2 proc~dodpc3 dodpc3 proc~dodpcr->proc~dodpc3 proc~dodphd dodphd proc~dodpcr->proc~dodphd proc~dhstep dhstep proc~dodpc1->proc~dhstep proc~dppt dppt proc~dodpc3->proc~dppt proc~dppnml dppnml proc~dppt->proc~dppnml

Called by

proc~~dodpcr~~CalledByGraph proc~dodpcr dodpcr proc~dodmn dodmn proc~dodmn->proc~dodpcr proc~doddrv doddrv proc~doddrv->proc~dodmn 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
real(kind=wp), public :: pnlty
logical, public :: anajac
logical, public :: cdjac
logical, public :: chkjac
logical, public :: dovcv
logical, public :: implct
logical, public :: initd
logical, public :: isodr
logical, public :: redoj
logical, public :: restrt
character(len=3), public :: typ

Source Code

   impure subroutine dodpcr &
      (ipr, lunrpt, &
      head, prtpen, fstitr, didvcv, iflag, &
      n, m, np, nq, npp, nnzw, &
      msgb, msgd, beta, y, ldy, x, ldx, delta, &
      we, ldwe, ld2we, wd, ldwd, ld2wd, &
      ifixb, ifixx, ldifx, &
      lower, upper, &
      ssf, tt, ldtt, stpb, stpd, ldstpd, &
      job, neta, taufac, sstol, partol, maxit, &
      wss, rvar, idf, sdbeta, &
      niter, nfev, njev, actred, prered, &
      tau, pnorm, alpha, f, rcond, irank, info, istop)
   !! Generate computation reports.
   ! Routines Called  DFLAGS, DODPC1, DODPC2, DODPC3, DODPHD
   ! Date Written   860529   (YYMMDD)
   ! Revision Date  920619   (YYMMDD)

      use odrpack_core, only: dflags

      integer, intent(in) :: ipr
         !! The variable indicating what is to be printed.
      integer, intent(in) :: lunrpt
         !! The logical unit number used for computation reports.
      logical, intent(inout) :: head
         !! The variable designating whether the heading is to be printed (`head = .true.`)
         !! or not (`head = .false.`).
      logical, intent(in) :: prtpen
         !! The variable designating whether the penalty parameter is to be printed in the
         !! iteration report (`prtpen = .true.`) or not (`prtpen = .false.`).
      logical, intent(in) :: fstitr
         !! The variable designating whether this is the first iteration (`fstitr = .true.`)
         !! or not (`fstitr = .false.`).
      logical, intent(in) :: didvcv
         !! The variable designating whether the covariance matrix was computed
         !! (`didvcv = .true.`) or not (`didvcv = .false.`).
      integer, intent(in) :: iflag
         !! The variable designating what is to be printed.
      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.
      integer, intent(in) :: npp
         !! The number of function parameters being estimated.
      integer, intent(in) :: nnzw
         !! The number of nonzero weighted observations.
      integer, intent(in) :: msgb(nq*np + 1)
         !! The error checking results for the Jacobian with respect to `beta`.
      integer, intent(in) :: msgd(nq*m + 1)
         !! The error checking results for the Jacobian with respect to `delta`.
      real(wp), intent(in) :: beta(np)
         !! The function parameters.
      real(wp), intent(in) :: y(ldy, nq)
         !! The response variable. Unused when the model is implicit.
      integer, intent(in) :: ldy
         !! The leading dimension of array `y`.
      real(wp), intent(in) :: x(ldx, m)
         !! The explanatory variable.
      integer, intent(in) :: ldx
         !! The leading dimension of array `x`.
      real(wp), intent(in) :: delta(n, m)
         !! The estimated errors in the explanatory variables.
      real(wp), intent(in) :: we(ldwe, ld2we, nq)
         !! The `epsilon` weights.
      integer, intent(in) :: ldwe
         !! The leading dimension of array `we`.
      integer, intent(in) :: ld2we
         !! The second dimension of array `we`.
      real(wp), intent(in) :: wd(ldwd, ld2wd, m)
         !! The `delta` weights.
      integer, intent(in) :: ldwd
         !! The leading dimension of array `wd`.
      integer, intent(in) :: ld2wd
         !! The second dimension of array `wd`.
      integer, intent(in) :: ifixb(np)
         !! The values designating whether the elements of `beta` are fixed at their input 
         !! values or not.
      integer, intent(in) :: ifixx(ldifx, m)
         !! The values designating whether the elements of `x` are fixed at their input values
         !! or not.
      integer, intent(in) :: ldifx
         !! The leading dimension of array `ifixx`.
      real(wp), intent(in) :: lower(np)
         !! The lower bounds for `beta`.
      real(wp), intent(in) :: upper(np)
         !! The upper bounds for `beta`.
      real(wp), intent(in) :: ssf(np)
         !! The scaling values for `beta`.
      real(wp), intent(in) :: tt(ldtt, m)
         !! The scaling values for `delta`.
      integer, intent(in) :: ldtt
         !! The leading dimension of array `tt`.
      real(wp), intent(in) :: stpb(np)
         !! The relative step used for computing finite difference derivatives with respect
         !! to `beta`.
      real(wp), intent(in) :: stpd(ldstpd, m)
         !! The relative step used for computing finite difference derivatives with respect
         !! to `delta`.
      integer, intent(in) :: ldstpd
         !! The leading dimension of array `stpd`.
      integer, intent(in) :: job
         !! The variable controlling problem initialization and computational method.
      integer, intent(in) :: neta
         !! The number of accurate digits in the function results.
      real(wp), intent(in) :: taufac
         !! The factor used to compute the initial trust region diameter.
      real(wp), intent(in) :: sstol
         !! The sum-of-squares convergence stopping tolerance.
      real(wp), intent(in) :: partol
         !! The parameter convergence stopping tolerance.
      integer, intent(in) :: maxit
         !! The maximum number of iterations allowed.
      real(wp), intent(in) :: wss(3)
         !! The sum-of-squares of the weighted `epsilon`s and `delta`s, the sum-of-squares of 
         !! the weighted `delta`s, and the sum-of-squares of the weighted `epsilon`s.
      real(wp), intent(in) :: rvar
         !! The residual variance.
      integer, intent(in) :: idf
         !! The degrees of freedom of the fit, equal to the number of observations with nonzero
         !! weighted derivatives minus the number of parameters being estimated.
      real(wp), intent(in) :: sdbeta(np)
         !! The standard errors of the estimated parameters.
      integer, intent(in) :: niter
         !! The number of iterations.
      integer, intent(in) :: nfev
         !! The number of function evaluations.
      integer, intent(in) :: njev
         !! The number of Jacobian evaluations.
      real(wp), intent(in) :: actred
         !! The actual relative reduction in the sum-of-squares.
      real(wp), intent(in) :: prered
         !! The predicted relative reduction in the sum-of-squares.
      real(wp), intent(in) :: tau
         !! The trust region diameter.
      real(wp), intent(in) :: pnorm
         !! The norm of the scaled estimated parameters.
      real(wp), intent(in) :: alpha
         !! The Levenberg-Marquardt parameter.
      real(wp), intent(in) :: f(n, nq)
         !! The estimated values of `epsilon`.
      real(wp), intent(in) :: rcond
         !! The approximate reciprocal condition of `tfjacb`.
      integer, intent(in) :: irank
         !! The rank deficiency of the Jacobian with respect to `beta`.
      integer, intent(in) :: info
         !! The variable designating why the computations were stopped.
      integer, intent(in) :: istop
         !! The variable designating whether there are problems computing the function at the
         !! current `beta` and `delta`.

      ! Local scalars
      real(wp) :: pnlty
      logical :: anajac, cdjac, chkjac, dovcv, implct, initd, isodr, redoj, restrt
      character(len=3) :: typ

      ! Variable Definitions (alphabetically)
      !  ACTRED:  The actual relative reduction in the sum-of-squares.
      !  ALPHA:   The Levenberg-Marquardt parameter.
      !  ANAJAC:  The variable designating whether the Jacobians are computed by finite
      !           differences (ANAJAC=FALSE) or not (ANAJAC=TRUE).
      !  BETA:    The function parameters.
      !  CDJAC:   The variable designating whether the jacobians are computed by central
      !           differences (CDJAC=TRUE) or by forward differences (CDJAC=FALSE).
      !  CHKJAC:  The variable designating whether the user supplied Jacobians are to be checked
      !           (CHKJAC=TRUE) or not (CHKJAC=FALSE).
      !  DELTA:   The estimated errors in the explanatory variables.
      !  DIDVCV:  The variable designating whether the covariance matrix was computed
      !           (DIDVCV=TRUE) or not (DIDVCV=FALSE).
      !  DOVCV:   The variable designating whether the covariance matrix was to be computed
      !           (DOVCV=TRUE) or not (DOVCV=FALSE).
      !  F:       The (weighted) estimated values of EPSILON.
      !  FSTITR:  The variable designating whether this is the first iteration (FSTITR=TRUE) or
      !           not (FSTITR=FALSE).
      !  HEAD:    The variable designating whether the heading is to be printed (HEAD=TRUE) or
      !           not (HEAD=FALSE).
      !  IDF:     The degrees of freedom of the fit, equal to the number of observations with
      !           nonzero weighted derivatives minus the number of parameters being estimated.
      !  IFIXB:   The values designating whether the elements of BETA are fixed at their input
      !           values or not.
      !  IFIXX:   The values designating whether the elements of X are fixed at their input
      !           values or not.
      !  IFLAG:   The variable designating what is to be printed.
      !  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.
      !  INITD:   The variable designating whether DELTA is initialized to zero (INITD=TRUE) or
      !           to the values in the first N by M elements of array WORK (INITD=FALSE).
      !  IPR:     The value indicating the report to be printed.
      !  IRANK:   The rank deficiency of the Jacobian wrt BETA.
      !  ISODR:   The variable designating whether the solution is by ODR (ISODR=TRUE) or by
      !           OLS (ISODR=FALSE).
      !  ISTOP:   The variable designating whether there are problems computing the function at
      !           the current BETA and DELTA.
      !  JOB:     The variable controling problem initialization and computational method.
      !  LDIFX:   The leading dimension of array IFIXX.
      !  LDSTPD:  The leading dimension of array STPD.
      !  LDTT:    The leading dimension of array TT.
      !  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 Y.
      !  LD2WD:   The second dimension of array WD.
      !  LD2WE:   The second dimension of array WE.
      !  LOWER:   Lower bound on BETA.
      !  LUNRPT:  The logical unit number for computation reports.
      !  M:       The number of columns of data in the explanatory variable.
      !  MAXIT:   The maximum number of iterations allowed.
      !  MSGB:    The error checking results for the Jacobian wrt BETA.
      !  MSGD:    The error checking results for the Jacobian wrt DELTA.
      !  N:       The number of observations.
      !  NETA:    The number of accurate digits in the function results.
      !  NFEV:    The number of function evaluations.
      !  NITER:   The number of iterations.
      !  NJEV:    The number of Jacobian evaluations.
      !  NNZW:    The number of nonzero weighted observations.
      !  NP:      The number of function parameters.
      !  NQ:      The number of responses per observation.
      !  NPP:     The number of function parameters being estimated.
      !  PARTOL:  The parameter convergence stopping tolerance.
      !  PNLTY:   The penalty parameter for an implicit model.
      !  PNORM:   The norm of the scaled estimated parameters.
      !  PRERED:  The predicted relative reduction in the sum-of-squares.
      !  PRTPEN:  The variable designating whether the penalty parameter is  to be printed in 
      !           the iteration report (PRTPEN=TRUE) or not (PRTPEN=FALSE).
      !  RCOND:   The approximate reciprocal condition number of TFJACB.
      !  REDOJ:   The variable designating whether the Jacobian matrix is to be recomputed for 
      !           the computation of the covariance matrix (REDOJ=TRUE) or not (REDOJ=FALSE).
      !  RESTRT:  The variable designating whether the call is a restart (RESTRT=TRUE) or not
      !           (RESTRT=FALSE).
      !  RVAR:    The residual variance.
      !  SDBETA:  The standard deviations of the estimated BETA'S.
      !  SSF:     The scaling values for BETA.
      !  SSTOL:   The sum-of-squares convergence stopping tolerance.
      !  STPB:    The relative step for computing finite difference derivatives with respect
      !           to BETA.
      !  STPD:    The relative step for computing finite difference derivatives with respect
      !           to DELTA.
      !  TAU:     The trust region diameter.
      !  TAUFAC:  The factor used to compute the initial trust region diameter.
      !  TT:      The scaling values for DELTA.
      !  TYP:     The CHARACTER*3 string "ODR" or "OLS".
      !  UPPER:   Upper bound on BETA.
      !  WE:      The EPSILON weights.
      !  WD:      The DELTA weights.
      !  WSS:     The sum-of-squares of the weighted EPSILONS and DELTAS, the sum-of-squares of
      !           the weighted DELTAS, and the sum-of-squares of the weighted EPSILONS.
      !  X:       The explanatory variable.
      !  Y:       The dependent variable.  Unused when the model is implicit.

      call dflags(job, restrt, initd, dovcv, redoj, anajac, cdjac, chkjac, isodr, implct)
      pnlty = abs(we(1, 1, 1))

      if (head) then
         call dodphd(head, lunrpt)
      end if
      if (isodr) then
         typ = 'ODR'
      else
         typ = 'OLS'
      end if

      ! Print initial summary
      if (iflag == 1) then
         write (lunrpt, 1200) typ
         call dodpc1 &
            (ipr, lunrpt, &
            anajac, cdjac, chkjac, initd, restrt, isodr, implct, dovcv, redoj, &
            msgb(1), msgb(2), msgd(1), msgd(2), n, m, np, nq, npp, nnzw, x, &
            ldx, ifixx, ldifx, delta, wd, ldwd, ld2wd, tt, ldtt, stpd, ldstpd, &
            y, ldy, we, ldwe, ld2we, pnlty, beta, ifixb, ssf, stpb, lower, &
            upper, job, neta, taufac, sstol, partol, maxit, wss(1), wss(2), &
            wss(3))

         ! Print iteration reports
      elseif (iflag == 2) then

         if (fstitr) then
            write (lunrpt, 1300) typ
         end if
         call dodpc2 &
            (ipr, lunrpt, fstitr, implct, prtpen, &
            pnlty, &
            niter, nfev, wss(1), actred, prered, alpha, tau, pnorm, np, beta)

         ! Print final summary
      elseif (iflag == 3) then

         write (lunrpt, 1400) typ
         call dodpc3 &
            (ipr, lunrpt, &
            isodr, implct, didvcv, dovcv, redoj, anajac, &
            n, m, np, nq, npp, &
            info, niter, nfev, njev, irank, rcond, istop, &
            wss(1), wss(2), wss(3), pnlty, rvar, idf, &
            beta, sdbeta, ifixb, f, delta, lower, upper)
      end if

      ! Format statements

   1200 format &
         (/' *** Initial summary for fit by method of ', A3, ' ***')
   1300 format &
         (/' *** Iteration reports for fit by method of ', A3, ' ***')
   1400 format &
         (/' *** Final summary for fit by method of ', A3, ' ***')

   end subroutine dodpcr