dodpc3 Subroutine

public impure subroutine dodpc3(ipr, lunrpt, isodr, implct, didvcv, dovcv, redoj, anajac, n, m, np, nq, npp, info, niter, nfev, njev, irank, rcond, istop, wss, wssdel, wsseps, pnlty, rvar, idf, beta, sdbeta, ifixb2, f, delta, lower, upper)

Uses

  • proc~~dodpc3~~UsesGraph proc~dodpc3 dodpc3 module~odrpack_core odrpack_core proc~dodpc3->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 final summary report.

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(in) :: isodr

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

logical, intent(in) :: implct

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

logical, intent(in) :: didvcv

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

logical, intent(in) :: dovcv

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

logical, intent(in) :: 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.).

logical, intent(in) :: anajac

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

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) :: info

The variable designating why the computations were stopped.

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.

integer, intent(in) :: irank

The rank deficiency of the Jacobian with respect to beta.

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

The approximate reciprocal condition of tfjacb.

integer, intent(in) :: istop

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

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

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

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

The sum-of-squares of the weighted deltas.

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

The sum-of-squares of the weighted epsilons.

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

The penalty parameter for an implicit model.

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) :: beta(np)

The function parameters.

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

The standard errors of the estimated parameters.

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

The values designating whether the elements of beta were estimated, fixed, or dropped because they caused rank deficiency.

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

The estimated values of epsilon.

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

The estimated errors in the explanatory variables.

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

Lower bound on beta.

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

Upper bound on beta.


Calls

proc~~dodpc3~~CallsGraph proc~dodpc3 dodpc3 proc~dppt dppt proc~dodpc3->proc~dppt proc~dppnml dppnml proc~dppt->proc~dppnml

Called by

proc~~dodpc3~~CalledByGraph proc~dodpc3 dodpc3 proc~dodpcr dodpcr proc~dodpcr->proc~dodpc3 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 :: tval
integer, public :: d1
integer, public :: d2
integer, public :: d3
integer, public :: d4
integer, public :: d5
integer, public :: i
integer, public :: j
integer, public :: k
integer, public :: l
integer, public :: nplm1
character(len=90), public :: fmt1

Source Code

   impure subroutine dodpc3 &
      (ipr, lunrpt, &
      isodr, implct, didvcv, dovcv, redoj, anajac, &
      n, m, np, nq, npp, &
      info, niter, nfev, njev, irank, rcond, istop, &
      wss, wssdel, wsseps, pnlty, rvar, idf, &
      beta, sdbeta, ifixb2, f, delta, &
      lower, upper)
   !! Generate final summary report.
   ! Routines Called  DPPT
   ! Date Written   860529   (YYMMDD)
   ! REvision Date  920619   (YYMMDD)

      use odrpack_core, only: dppt

      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(in) :: isodr
         !! The variable designating whether the solution is by ODR (`isodr = .true.`) or
         !! by OLS (`isodr = .false.`).
      logical, intent(in) :: implct
         !! The variable designating whether the solution is by implicit ODR (`implct = .true.`)
         !! or explicit ODR (`implct = .false.`).
      logical, intent(in) :: didvcv
         !! The variable designating whether the covariance matrix was computed (`didvcv = .true.`)
         !! or not (`didvcv = .false.`).
      logical, intent(in) :: dovcv
         !! The variable designating whether the covariance matrix was to be computed
         !! (`dovcv = .true.`) or not (`dovcv = .false.`).
      logical, intent(in) :: 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.`).
      logical, intent(in) :: anajac
         !! The variable designating whether the Jacobians are computed by finite differences
         !! (`anajac = .false.`) or not (`anajac = .true.`).
      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) :: info
         !! The variable designating why the computations were stopped.
      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.
      integer, intent(in) :: irank
         !! The rank deficiency of the Jacobian with respect to `beta`.
      real(wp), intent(in) :: rcond
         !! The approximate reciprocal condition of `tfjacb`.
      integer, intent(in) :: istop
         !! The variable designating whether there are problems computing the function at the
         !! current `beta` and `delta`.
      real(wp), intent(in) :: wss
         !! The sum-of-squares of the weighted `epsilon`s and `delta`s.
      real(wp), intent(in) :: wssdel
         !! The sum-of-squares of the weighted `delta`s.
      real(wp), intent(in) :: wsseps
         !! The sum-of-squares of the weighted `epsilon`s.
      real(wp), intent(in) :: pnlty
         !! The penalty parameter for an implicit model.
      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) :: beta(np)
         !! The function parameters.
      real(wp), intent(in) :: sdbeta(np)
         !! The standard errors of the estimated parameters.
      integer, intent(in) :: ifixb2(np)
         !! The values designating whether the elements of `beta` were estimated, fixed, or
         !! dropped because they caused rank deficiency.
      real(wp), intent(in) :: f(n, nq)
         !! The estimated values of `epsilon`.
      real(wp), intent(in) :: delta(n, m)
         !! The estimated errors in the explanatory variables.
      real(wp), intent(in) :: lower(np)
         !! Lower bound on `beta`.
      real(wp), intent(in) :: upper(np)
         !! Upper bound on `beta`.

      ! Local scalars
      real(wp) :: tval
      integer :: d1, d2, d3, d4, d5, i, j, k, l, nplm1
      character(len=90) :: fmt1

      ! Variable Definitions (alphabetically)
      !  ANAJAC:  The variable designating whether the JACOBIANS are computed by finite
      !           differences (ANAJAC=FALSE) or not (ANAJAC=TRUE).
      !  BETA:    The function parameters.
      !  D1:      The first digit of INFO.
      !  D2:      The second digit of INFO.
      !  D3:      The third digit of INFO.
      !  D4:      The fourth digit of INFO.
      !  D5:      The fifth digit of INFO.
      !  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 estimated values of EPSILON.
      !  FMT1:    A CHARACTER*90 variable used for formats.
      !  I:       An indexing variable.
      !  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.
      !  IFIXB2:  The values designating whether the elements of BETA were estimated, fixed,
      !           or dropped because they caused rank deficiency, corresponding to values of
      !           IFIXB2 equaling 1, 0, and -1, respectively.  If IFIXB2 is -2, then no attempt
      !           was made to estimate the parameters because MAXIT = 0.
      !  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.
      !  IPR:     The variable indicating what is 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.
      !  J:       An indexing variable.
      !  K:       An indexing variable.
      !  L:       An indexing variable.
      !  LOWER:   Lower bound on BETA.
      !  LUNRPT:  The logical unit number used for computation reports.
      !  M:       The number of columns of data in the explanatory variable.
      !  N:       The number of observations.
      !  NFEV:    The number of function evaluations.
      !  NITER:   The number of iterations.
      !  NJEV:    The number of Jacobian evaluations.
      !  NP:      The number of function parameters.
      !  NPLM1:   The number of items to be printed per line, minus one.
      !  NPP:     The number of function parameters being estimated.
      !  NQ:      The number of responses per observation.
      !  PNLTY:   The penalty parameter for an implicit model.
      !  RCOND:   The approximate reciprocal condition 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).
      !  RVAR:    The residual variance.
      !  SDBETA:  The standard errors of the estimated parameters.
      !  TVAL:    The value of the 97.5 percent point function for the T distribution.
      !  UPPER:   Upper bound on BETA.
      !  WSS:     The sum-of-squares of the weighted EPSILONS and DELTAS.
      !  WSSDEL:  The sum-of-squares of the weighted DELTAS.
      !  WSSEPS:  The sum-of-squares of the weighted EPSILONS.

      d1 = info/10000
      d2 = mod(info, 10000)/1000
      d3 = mod(info, 1000)/100
      d4 = mod(info, 100)/10
      d5 = mod(info, 10)

      ! Print stopping conditions
      write (lunrpt, 1000)
      if (info <= 9) then
         if (info == 1) then
            write (lunrpt, 1011) info
         elseif (info == 2) then
            write (lunrpt, 1012) info
         elseif (info == 3) then
            write (lunrpt, 1013) info
         elseif (info == 4) then
            write (lunrpt, 1014) info
         elseif (info <= 9) then
            write (lunrpt, 1015) info
         end if
      elseif (info <= 9999) then

         ! Print warning diagnostics
         write (lunrpt, 1020) info
         if (d2 == 1) write (lunrpt, 1021)
         if (d3 == 1) write (lunrpt, 1022)
         if (d4 == 1) write (lunrpt, 1023)
         if (d4 == 2) write (lunrpt, 1024)
         if (d5 == 1) then
            write (lunrpt, 1031)
         elseif (d5 == 2) then
            write (lunrpt, 1032)
         elseif (d5 == 3) then
            write (lunrpt, 1033)
         elseif (d5 == 4) then
            write (lunrpt, 1034)
         elseif (d5 <= 9) then
            write (lunrpt, 1035) d5
         end if
      else

         ! Print error messages
         write (lunrpt, 1040) info
         if (d1 == 5) then
            write (lunrpt, 1042)
            if (d2 /= 0) write (lunrpt, 1043) d2
            if (d3 == 3) then
               write (lunrpt, 1044) d3
            elseif (d3 /= 0) then
               write (lunrpt, 1045) d3
            end if
         elseif (d1 == 6) then
            write (lunrpt, 1050)
         else
            write (lunrpt, 1060) d1
         end if
      end if

      ! Print misc. stopping info
      write (lunrpt, 1300) niter
      write (lunrpt, 1310) nfev
      if (anajac) write (lunrpt, 1320) njev
      write (lunrpt, 1330) irank
      write (lunrpt, 1340) rcond
      write (lunrpt, 1350) istop

      ! Print final sum of squares
      if (implct) then
         write (lunrpt, 2000) wssdel
         if (isodr) then
            write (lunrpt, 2010) wss, wsseps, pnlty
         end if
      else
         write (lunrpt, 2100) wss
         if (isodr) then
            write (lunrpt, 2110) wssdel, wsseps
         end if
      end if
      if (didvcv) then
         write (lunrpt, 2200) sqrt(rvar), idf
      end if

      nplm1 = 3

      ! Print estimated BETA's, and, if, full rank, their standard errors
      write (lunrpt, 3000)
      if (didvcv) then
         write (lunrpt, 7300)
         tval = dppt(0.975E0_wp, idf)
         do j = 1, np
            if (ifixb2(j) >= 1) then
               write (lunrpt, 8400) j, beta(j), &
                  lower(j), upper(j), &
                  sdbeta(j), &
                  beta(j) - tval*sdbeta(j), &
                  beta(j) + tval*sdbeta(j)
            elseif (ifixb2(j) == 0) then
               write (lunrpt, 8600) j, beta(j), lower(j), upper(j)
            else
               write (lunrpt, 8700) j, beta(j), lower(j), upper(j)
            end if
         end do
         if (.not. redoj) write (lunrpt, 7310)
      else
         if (dovcv) then
            if (d1 <= 5) then
               write (lunrpt, 7410)
            else
               write (lunrpt, 7420)
            end if
         end if

         if ((irank == 0 .and. npp == np) .or. niter == 0) then
            if (np == 1) then
               write (lunrpt, 7100)
            else
               write (lunrpt, 7200)
            end if
            do j = 1, np, nplm1 + 1
               k = min(j + nplm1, np)
               if (k == j) then
                  write (lunrpt, 8100) j, beta(j)
               else
                  write (lunrpt, 8200) j, k, (beta(l), l=j, k)
               end if
            end do
            if (niter >= 1) then
               write (lunrpt, 8800)
            else
               write (lunrpt, 8900)
            end if
         else
            write (lunrpt, 7500)
            do j = 1, np
               if (ifixb2(j) >= 1) then
                  write (lunrpt, 8500) j, beta(j), lower(j), upper(j)
               elseif (ifixb2(j) == 0) then
                  write (lunrpt, 8600) j, beta(j), lower(j), upper(j)
               else
                  write (lunrpt, 8700) j, beta(j), lower(j), upper(j)
               end if
            end do
         end if
      end if

      if (ipr == 1) return

      ! Print EPSILON's and DELTA's together in a column if the number of
      ! columns of data in EPSILON and DELTA is less than or equal to three.
      if (implct .and. &
         (m <= 4)) then
         write (lunrpt, 4100)
         write (fmt1, 9110) m
         write (lunrpt, fmt1) (j, j=1, m)
         do i = 1, n
            write (lunrpt, 4130) i, (delta(i, j), j=1, m)
         end do

      elseif (isodr .and. &
            (nq + m <= 4)) then
         write (lunrpt, 4110)
         write (fmt1, 9120) nq, m
         write (lunrpt, fmt1) (l, l=1, nq), (j, j=1, m)
         do i = 1, n
            write (lunrpt, 4130) i, (f(i, l), l=1, nq), (delta(i, j), j=1, m)
         end do

      elseif (.not. isodr .and. &
            ((nq >= 2) .and. &
               (nq <= 4))) then
         write (lunrpt, 4120)
         write (fmt1, 9130) nq
         write (lunrpt, fmt1) (l, l=1, nq)
         do i = 1, n
            write (lunrpt, 4130) i, (f(i, l), l=1, nq)
         end do
      else

         ! Print EPSILON's and DELTA's separately
         if (.not. implct) then

            ! Print EPSILON'S
            do j = 1, nq
               write (lunrpt, 4200) j
               if (n == 1) then
                  write (lunrpt, 7100)
               else
                  write (lunrpt, 7200)
               end if
               do i = 1, n, nplm1 + 1
                  k = min(i + nplm1, n)
                  if (i == k) then
                     write (lunrpt, 8100) i, f(i, j)
                  else
                     write (lunrpt, 8200) i, k, (f(l, j), l=i, k)
                  end if
               end do
            end do
         end if

         ! Print DELTA'S
         if (isodr) then
            do j = 1, m
               write (lunrpt, 4300) j
               if (n == 1) then
                  write (lunrpt, 7100)
               else
                  write (lunrpt, 7200)
               end if
               do i = 1, n, nplm1 + 1
                  k = min(i + nplm1, n)
                  if (i == k) then
                     write (lunrpt, 8100) i, delta(i, j)
                  else
                     write (lunrpt, 8200) i, k, (delta(l, j), l=i, k)
                  end if
               end do
            end do
         end if
      end if

      ! Format statements

   1000 format &
         (/' --- Stopping Conditions:')
   1011 format &
         ('         INFO = ', I5, ' ==> sum of squares convergence.')
   1012 format &
         ('         INFO = ', I5, ' ==> parameter convergence.')
   1013 format &
         ('         INFO = ', I5, ' ==> sum of squares convergence and', &
         ' parameter convergence.')
   1014 format &
         ('         INFO = ', I5, ' ==> iteration limit reached.')
   1015 format &
         ('         INFO = ', I5, ' ==> unexpected value,', &
         ' probably indicating'/ &
         '                           incorrectly specified', &
         ' user input.')
   1020 format &
         ('         INFO = ', I5.4/ &
         '              =  ABCD, where a nonzero value for digit A,', &
         ' B, or C indicates why'/ &
         '                       the results might be questionable,', &
         ' and digit D indicates'/ &
         '                       the actual stopping condition.')
   1021 format &
         ('                       A=1 ==> derivatives are', &
         ' questionable.')
   1022 format &
         ('                       B=1 ==> user set ISTOP to', &
         ' nonzero value during last'/ &
         '                               call to subroutine FCN.')
   1023 format &
         ('                       C=1 ==> derivatives are not', &
         ' full rank at the solution.')
   1024 format &
         ('                       C=2 ==> derivatives are zero', &
         ' rank at the solution.')
   1031 format &
         ('                       D=1 ==> sum of squares convergence.')
   1032 format &
         ('                       D=2 ==> parameter convergence.')
   1033 format &
         ('                       D=3 ==> sum of squares convergence', &
         ' and parameter convergence.')
   1034 format &
         ('                       D=4 ==> iteration limit reached.')
   1035 format &
         ('                       D=', I1, ' ==> unexpected value,', &
         ' probably indicating'/ &
         '                               incorrectly specified', &
         ' user input.')
   1040 format &
         ('         INFO = ', I5.5/ &
         '              = ABCDE, where a nonzero value for a given', &
         ' digit indicates an'/ &
         '                       abnormal stopping condition.')
   1042 format &
         ('                       A=5 ==> user stopped computations', &
         ' in subroutine FCN.')
   1043 format &
         ('                       B=', I1, ' ==> computations were', &
         ' stopped during the'/ &
         '                                    function evaluation.')
   1044 format &
         ('                       C=', I1, ' ==> computations were', &
         ' stopped because'/ &
         '                                    derivatives with', &
         ' respect to delta were'/ &
         '                                    computed by', &
         ' subroutine FCN when'/ &
         '                                    fit is OLS.')
   1045 format &
         ('                       C=', I1, ' ==> computations were', &
         ' stopped during the'/ &
         '                                    jacobian evaluation.')
   1050 format &
         ('                       A=6 ==> numerical instabilities', &
         ' have been detected,'/ &
         '                               possibly indicating', &
         ' a discontinuity in the'/ &
         '                               derivatives or a poor', &
         ' poor choice of problem'/ &
         '                               scale or weights.')
   1060 format &
         ('                       A=', I1, ' ==> unexpected value,', &
         ' probably indicating'/ &
         '                               incorrectly specified', &
         ' user input.')
   1300 format &
         ('        NITER = ', I5, &
         '          (number of iterations)')
   1310 format &
         ('         NFEV = ', I5, &
         '          (number of function evaluations)')
   1320 format &
         ('         NJEV = ', I5, &
         '          (number of jacobian evaluations)')
   1330 format &
         ('        IRANK = ', I5, &
         '          (rank deficiency)')
   1340 format &
         ('        RCOND = ', 1P, E12.2, &
         '   (inverse condition number)')
   !1341 FORMAT
   !       +  ('                      ==> POSSIBLY FEWER THAN 2 SIGNIFICANT',
   !       +                        ' DIGITS IN RESULTS;'/
   !       +   '                          SEE ODRPACK95 REFERENCE',
   !       +                        ' GUIDE, SECTION 4.C.')
   1350 format &
         ('        ISTOP = ', I5, &
         '          (returned by user from', &
         ' subroutine FCN)')
   2000 format &
         (/' --- Final Sum of Squared Weighted Deltas = ', &
         17X, 1P, E17.8)
   2010 format &
         ('         Final Penalty Function Value     = ', 1P, E17.8/ &
         '               Penalty Term               = ', 1P, E17.8/ &
         '               Penalty Parameter          = ', 1P, E10.1)
   2100 format &
         (/' --- Final Weighted Sums of Squares       = ', 17X, 1P, E17.8)
   2110 format &
         ('         Sum of Squared Weighted Deltas   = ', 1P, E17.8/ &
         '         Sum of Squared Weighted Epsilons = ', 1P, E17.8)
   2200 format &
         (/' --- Residual Standard Deviation          = ', &
         17X, 1P, E17.8/ &
         '         Degrees of Freedom               =', I5)
   3000 format &
         (/' --- Estimated BETA(J), J = 1, ..., NP:')
   4100 format &
         (/' --- Estimated DELTA(I,*), I = 1, ..., N:')
   4110 format &
         (/' --- Estimated EPSILON(I) and DELTA(I,*), I = 1, ..., N:')
   4120 format &
         (/' --- Estimated EPSILON(I), I = 1, ..., N:')
   4130 format(5X, I5, 1P, 5E16.8)
   4200 format &
         (/' --- Estimated EPSILON(I,', I3, '), I = 1, ..., N:')
   4300 format &
         (/' --- Estimated DELTA(I,', I3, '), I = 1, ..., N:')
   7100 format &
         (/'           Index           Value'/)
   7200 format &
         (/'           Index           Value -------------->'/)
   7300 format &
         (/'                     BETA      LOWER     UPPER      S.D. ', &
         ' ___ 95% Confidence ___'/ &
         '                                                    BETA ', &
         '        Interval'/)
   7310 format &
         (/'     N.B. standard errors and confidence intervals are', &
         ' computed using'/ &
         '          derivatives calculated at the beginning', &
         ' of the last iteration,'/ &
         '          and not using derivatives re-evaluated at the', &
         ' final solution.')
   7410 format &
         (/'     N.B. the standard errors of the estimated betas were', &
         ' not computed because'/ &
         '          the derivatives were not available.  Either MAXIT', &
         ' is 0 and the third'/ &
         '          digit of JOB is greater than 1, or the most', &
         ' recently tried values of'/ &
         '          BETA and/or X+DELTA were identified as', &
         ' unacceptable by user supplied'/ &
         '          subroutine FCN.')
   7420 format &
         (/'     N.B. the standard errors of the estimated betas were', &
         ' not computed.'/ &
         '          (see info above.)')
   7500 format &
         (/'                     BETA         Status')
   8100 format &
         (11X, I5, 1P, E16.8)
   8200 format &
         (3X, I5, ' to', I5, 1P, 7E16.8)
   8400 format &
         (3X, I5, 1X, 1P, E16.8, 1X, E10.2, E10.2, E10.2, 1X, E10.2, 1X, 'to', E10.2)
   8500 format &
         (3X, I5, 1X, 1P, E16.8, 1X, E10.2, E10.2, 4X, 'Estimated')
   8600 format &
         (3X, I5, 1X, 1P, E16.8, 1X, E10.2, E10.2, 4X, '    Fixed')
   8700 format &
         (3X, I5, 1X, 1P, E16.8, 1X, E10.2, E10.2, 4X, '  Dropped')
   8800 format &
         (/'     N.B. no parameters were fixed by the user or', &
         ' dropped at the last'/ &
         '          iteration because they caused the model to be', &
         ' rank deficient.')
   8900 format &
         (/'     N.B. no change was made to the user supplied parameter', &
         ' values because'/ &
         '          MAXIT=0.')
   9110 format &
         ('(/''         I'',', &
         I2, '(''      DELTA(I,'',I1,'')'')/)')
   9120 format &
         ('(/''         I'',', &
         I2, '(''    EPSILON(I,'',I1,'')''),', &
         I2, '(''      DELTA(I,'',I1,'')'')/)')
   9130 format &
         ('(/''         I'',', &
         I2, '(''    EPSILON(I,'',I1,'')'')/)')

   end subroutine dodpc3