print_reports Subroutine

public impure subroutine print_reports(ipr, lunrpt, head, printpen, firstitr, didvcv, iflag, n, m, np, q, npp, nnzw, msgb, msgd, beta, y, x, 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~~print_reports~~UsesGraph proc~print_reports print_reports module~odrpack_core odrpack_core proc~print_reports->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

Variable indicating what is to be printed.

integer, intent(in) :: lunrpt

Logical unit number used for computation reports.

logical, intent(inout) :: head

Variable designating whether the heading is to be printed (.true.) or not (.false.).

logical, intent(in) :: printpen

Variable designating whether the penalty parameter is to be printed in the iteration report (.true.) or not (.false.).

logical, intent(in) :: firstitr

Variable designating whether this is the first iteration (.true.) or not (.false.).

logical, intent(in) :: didvcv

Variable designating whether the covariance matrix was computed (.true.) or not (.false.).

integer, intent(in) :: iflag

Variable designating what is to be printed.

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.

integer, intent(in) :: npp

Number of function parameters being estimated.

integer, intent(in) :: nnzw

Number of nonzero weighted observations.

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

Error checking results for the Jacobian with respect to beta.

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

Error checking results for the Jacobian with respect to delta.

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

Function parameters.

real(kind=wp), intent(in) :: y(n,q)

Response variable. Unused when the model is implicit.

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

Explanatory variable.

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

Estimated errors in the explanatory variables.

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

epsilon weights.

integer, intent(in) :: ldwe

Leading dimension of array we.

integer, intent(in) :: ld2we

Second dimension of array we.

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

delta weights.

integer, intent(in) :: ldwd

Leading dimension of array wd.

integer, intent(in) :: ld2wd

Second dimension of array wd.

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

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

Lower bounds for beta.

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

Upper bounds for beta.

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

Scaling values for beta.

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

Scaling values for delta.

integer, intent(in) :: ldtt

Leading dimension of array tt.

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

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

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

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

integer, intent(in) :: ldstpd

Leading dimension of array stpd.

integer, intent(in) :: job

Variable controlling problem initialization and computational method.

integer, intent(in) :: neta

Number of accurate digits in the function results.

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

Factor used to compute the initial trust region diameter.

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

Sum-of-squares convergence stopping tolerance.

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

Parameter convergence stopping tolerance.

integer, intent(in) :: maxit

Maximum number of iterations allowed.

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

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

Residual variance.

integer, intent(in) :: idf

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)

Standard errors of the estimated parameters.

integer, intent(in) :: niter

Number of iterations.

integer, intent(in) :: nfev

Number of function evaluations.

integer, intent(in) :: njev

Number of Jacobian evaluations.

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

Actual relative reduction in the sum-of-squares.

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

Predicted relative reduction in the sum-of-squares.

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

Trust region diameter.

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

Norm of the scaled estimated parameters.

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

Levenberg-Marquardt parameter.

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

Estimated values of epsilon.

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

Approximate reciprocal condition of tfjacb.

integer, intent(in) :: irank

Rank deficiency of the Jacobian with respect to beta.

integer, intent(in) :: info

Variable designating why the computations were stopped.

integer, intent(in) :: istop

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


Calls

proc~~print_reports~~CallsGraph proc~print_reports print_reports proc~print_header print_header proc~print_reports->proc~print_header proc~print_report_final print_report_final proc~print_reports->proc~print_report_final proc~print_report_initial print_report_initial proc~print_reports->proc~print_report_initial proc~print_report_iter print_report_iter proc~print_reports->proc~print_report_iter proc~set_flags set_flags proc~print_reports->proc~set_flags proc~ppf_tstudent ppf_tstudent proc~print_report_final->proc~ppf_tstudent proc~hstep hstep proc~print_report_initial->proc~hstep proc~ppf_normal ppf_normal proc~ppf_tstudent->proc~ppf_normal

Called by

proc~~print_reports~~CalledByGraph proc~print_reports print_reports proc~odr odr proc~odr->proc~print_reports 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 :: penalty
logical, public :: anajac
logical, public :: cdjac
logical, public :: chkjac
logical, public :: dovcv
logical, public :: implicit
logical, public :: initd
logical, public :: isodr
logical, public :: redoj
logical, public :: restart
character(len=3), public :: typ

Source Code

   impure subroutine print_reports &
      (ipr, lunrpt, &
       head, printpen, firstitr, didvcv, iflag, &
       n, m, np, q, npp, nnzw, &
       msgb, msgd, beta, y, x, 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.

      use odrpack_core, only: set_flags

      integer, intent(in) :: ipr
         !! Variable indicating what is to be printed.
      integer, intent(in) :: lunrpt
         !! Logical unit number used for computation reports.
      logical, intent(inout) :: head
         !! Variable designating whether the heading is to be printed (`.true.`) or not (`.false.`).
      logical, intent(in) :: printpen
         !! Variable designating whether the penalty parameter is to be printed in the
         !! iteration report (`.true.`) or not (`.false.`).
      logical, intent(in) :: firstitr
         !! Variable designating whether this is the first iteration (`.true.`) or not (`.false.`).
      logical, intent(in) :: didvcv
         !! Variable designating whether the covariance matrix was computed (`.true.`)
         !! or not (`.false.`).
      integer, intent(in) :: iflag
         !! Variable designating what is to be printed.
      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.
      integer, intent(in) :: npp
         !! Number of function parameters being estimated.
      integer, intent(in) :: nnzw
         !! Number of nonzero weighted observations.
      integer, intent(in) :: msgb(q*np + 1)
         !! Error checking results for the Jacobian with respect to `beta`.
      integer, intent(in) :: msgd(q*m + 1)
         !! Error checking results for the Jacobian with respect to `delta`.
      real(wp), intent(in) :: beta(np)
         !! Function parameters.
      real(wp), intent(in) :: y(n, q)
         !! Response variable. Unused when the model is implicit.
      real(wp), intent(in) :: x(n, m)
         !! Explanatory variable.
      real(wp), intent(in) :: delta(n, m)
         !! Estimated errors in the explanatory variables.
      real(wp), intent(in) :: we(ldwe, ld2we, q)
         !! `epsilon` weights.
      integer, intent(in) :: ldwe
         !! Leading dimension of array `we`.
      integer, intent(in) :: ld2we
         !! Second dimension of array `we`.
      real(wp), intent(in) :: wd(ldwd, ld2wd, m)
         !! `delta` weights.
      integer, intent(in) :: ldwd
         !! Leading dimension of array `wd`.
      integer, intent(in) :: ld2wd
         !! Second dimension of array `wd`.
      integer, intent(in) :: ifixb(np)
         !! Values designating whether the elements of `beta` are fixed at their input
         !! values or not.
      integer, intent(in) :: ifixx(ldifx, m)
         !! Values designating whether the elements of `x` are fixed at their input values
         !! or not.
      integer, intent(in) :: ldifx
         !! Leading dimension of array `ifixx`.
      real(wp), intent(in) :: lower(np)
         !! Lower bounds for `beta`.
      real(wp), intent(in) :: upper(np)
         !! Upper bounds for `beta`.
      real(wp), intent(in) :: ssf(np)
         !! Scaling values for `beta`.
      real(wp), intent(in) :: tt(ldtt, m)
         !! Scaling values for `delta`.
      integer, intent(in) :: ldtt
         !! Leading dimension of array `tt`.
      real(wp), intent(in) :: stpb(np)
         !! Relative step used for computing finite difference derivatives with respect
         !! to `beta`.
      real(wp), intent(in) :: stpd(ldstpd, m)
         !! Relative step used for computing finite difference derivatives with respect
         !! to `delta`.
      integer, intent(in) :: ldstpd
         !! Leading dimension of array `stpd`.
      integer, intent(in) :: job
         !! Variable controlling problem initialization and computational method.
      integer, intent(in) :: neta
         !! Number of accurate digits in the function results.
      real(wp), intent(in) :: taufac
         !! Factor used to compute the initial trust region diameter.
      real(wp), intent(in) :: sstol
         !! Sum-of-squares convergence stopping tolerance.
      real(wp), intent(in) :: partol
         !! Parameter convergence stopping tolerance.
      integer, intent(in) :: maxit
         !! Maximum number of iterations allowed.
      real(wp), intent(in) :: wss(3)
         !! 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
         !! Residual variance.
      integer, intent(in) :: idf
         !! 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)
         !! Standard errors of the estimated parameters.
      integer, intent(in) :: niter
         !! Number of iterations.
      integer, intent(in) :: nfev
         !! Number of function evaluations.
      integer, intent(in) :: njev
         !! Number of Jacobian evaluations.
      real(wp), intent(in) :: actred
         !! Actual relative reduction in the sum-of-squares.
      real(wp), intent(in) :: prered
         !! Predicted relative reduction in the sum-of-squares.
      real(wp), intent(in) :: tau
         !! Trust region diameter.
      real(wp), intent(in) :: pnorm
         !! Norm of the scaled estimated parameters.
      real(wp), intent(in) :: alpha
         !! Levenberg-Marquardt parameter.
      real(wp), intent(in) :: f(n, q)
         !! Estimated values of `epsilon`.
      real(wp), intent(in) :: rcond
         !! Approximate reciprocal condition of `tfjacb`.
      integer, intent(in) :: irank
         !! Rank deficiency of the Jacobian with respect to `beta`.
      integer, intent(in) :: info
         !! Variable designating why the computations were stopped.
      integer, intent(in) :: istop
         !! Variable designating whether there are problems computing the function at the
         !! current `beta` and `delta`.

      ! Local scalars
      real(wp) :: penalty
      logical :: anajac, cdjac, chkjac, dovcv, implicit, initd, isodr, redoj, restart
      character(len=3) :: typ

      ! Variable Definitions (alphabetically)
      !  ANAJAC:   The variable designating whether the Jacobians are computed by finite
      !            differences (FALSE) or not (TRUE).
      !  CDJAC:    The variable designating whether the jacobians are computed by central
      !            differences (TRUE) or by forward differences (FALSE).
      !  CHKJAC:   The variable designating whether the user supplied Jacobians are to be checked
      !            (TRUE) or not (FALSE).
      !  DOVCV:    The variable designating whether the covariance matrix was to be computed
      !            (TRUE) or not (FALSE).
      !  IMPLICIT: The variable designating whether the solution is by implicit ODR (TRUE)
      !            or explicit ODR (FALSE).
      !  INITD:    The variable designating whether DELTA is initialized to zero (TRUE) or
      !            to the values in the first N by M elements of array RWORK (FALSE).
      !  ISODR:    The variable designating whether the solution is by ODR (TRUE) or by
      !            OLS (FALSE).
      !  PENALTY:  The penalty parameter for an implicit model.
      !  REDOJ:    The variable designating whether the Jacobian matrix is to be recomputed for
      !            the computation of the covariance matrix (TRUE) or not (FALSE).
      !  RESTART:  The variable designating whether the call is a restart (TRUE) or not
      !            (FALSE).
      !  TYP:      The CHARACTER*3 string "ODR" or "OLS".

      call set_flags(job, restart, initd, dovcv, redoj, anajac, cdjac, chkjac, isodr, implicit)
      penalty = abs(we(1, 1, 1))

      if (head) then
         call print_header(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 print_report_initial &
            (ipr, lunrpt, &
             anajac, cdjac, chkjac, initd, restart, isodr, implicit, dovcv, redoj, &
             msgb(1), msgb(2), msgd(1), msgd(2), n, m, np, q, npp, nnzw, x, &
             ifixx, ldifx, delta, wd, ldwd, ld2wd, tt, ldtt, stpd, ldstpd, &
             y, we, ldwe, ld2we, penalty, 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 (firstitr) then
            write (lunrpt, 1300) typ
         end if
         call print_report_iter &
            (ipr, lunrpt, firstitr, implicit, printpen, &
             penalty, &
             niter, nfev, wss(1), actred, prered, alpha, tau, pnorm, np, beta)

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

         write (lunrpt, 1400) typ
         call print_report_final &
            (ipr, lunrpt, &
             isodr, implicit, didvcv, dovcv, redoj, anajac, &
             n, m, np, q, npp, &
             info, niter, nfev, njev, irank, rcond, istop, &
             wss(1), wss(2), wss(3), penalty, 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 print_reports