odrpack_reports Module

Reporting routines.


Uses

  • module~~odrpack_reports~~UsesGraph module~odrpack_reports odrpack_reports module~odrpack_kinds odrpack_kinds module~odrpack_reports->module~odrpack_kinds iso_fortran_env iso_fortran_env module~odrpack_kinds->iso_fortran_env

Used by

  • module~~odrpack_reports~~UsedByGraph module~odrpack_reports odrpack_reports proc~odr odr proc~odr->module~odrpack_reports

Subroutines

public impure subroutine print_header(head, lunit)

Print report heading.

Arguments

Type IntentOptional Attributes Name
logical, intent(inout) :: head

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

integer, intent(in) :: lunit

Logical unit number to which the heading is written.

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)

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.

public impure subroutine print_report_initial(ipr, lunrpt, anajac, cdjac, chkjac, initd, restart, isodr, implicit, dovcv, redoj, msgb1, msgb, msgd1, msgd, 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, wssdel, wsseps)

Generate initial summary report.

Arguments

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

Value indicating the report to be printed.

integer, intent(in) :: lunrpt

Logical unit number for the computation reports.

logical, intent(in) :: anajac

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

logical, intent(in) :: cdjac

Variable designating whether the Jacobians are computed by central differences (.true.) or forward differences (.false.).

logical, intent(in) :: chkjac

Variable designating whether the user-supplied Jacobians are to be checked (.true.) or not (.false.).

logical, intent(in) :: initd

Variable designating whether delta is initialized to zero (.true.) or to the values in the first n by m elements of array rwork (.false.).

logical, intent(in) :: restart

Variable designating whether the call is a restart (.true.) or not (.false.).

logical, intent(in) :: isodr

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

logical, intent(in) :: implicit

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

logical, intent(in) :: dovcv

Variable designating whether the covariance matrix is to be computed (.true.) or not (.false.).

logical, intent(in) :: redoj

Variable designating whether the Jacobian matrix is to be recomputed for the computation of the covariance matrix (.true.) or not (.false.).

integer, intent(in) :: msgb1

Error checking results for the Jacobian with respect to beta.

integer, intent(in) :: msgb(q,np)

Error checking results for the Jacobian with respect to beta.

integer, intent(in) :: msgd1

Error checking results for the Jacobian with respect to delta.

integer, intent(in) :: msgd(q,m)

Error checking results for the Jacobian with respect to delta.

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 observational error weights.

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

Explanatory variable.

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) :: delta(n,m)

Estimated errors in the explanatory variables.

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.

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) :: stpd(ldstpd,m)

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

integer, intent(in) :: ldstpd

Leading dimension of array stpd.

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

Response variable. Unused when the model is implicit.

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

Penalty parameter for an implicit model.

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.

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

Scaling values for beta.

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

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

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

Lower bounds for beta.

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

Upper bounds for beta.

integer, intent(in) :: job

Variable controlling problem initialization and computational method.

integer, intent(in) :: neta

Number of accurate digits in the function results. A negative value indicates that neta was estimated by 'odrpack'. A positive value indicates the value was supplied by the user.

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

Sum-of-squares of the weighted epsilons and deltas.

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

Sum-of-squares of the weighted deltas.

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

Sum-of-squares of the weighted epsilons.

public impure subroutine print_report_iter(ipr, lunrpt, firstitr, implicit, printpen, penalty, niter, nfev, wss, actred, prered, alpha, tau, pnorm, np, beta)

Generate iteration reports.

Arguments

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

Value indicating the report to be printed.

integer, intent(in) :: lunrpt

Logical unit number used for computation reports.

logical, intent(in) :: firstitr

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

logical, intent(in) :: implicit

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

logical, intent(in) :: printpen

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

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

Penalty parameter for an implicit model.

integer, intent(in) :: niter

Number of iterations.

integer, intent(in) :: nfev

Number of function evaluations.

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

Sum-of-squares of the weighted epsilons and deltas.

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

Levenberg-Marquardt parameter.

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

Trust region diameter.

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

Norm of the scaled estimated parameters.

integer, intent(in) :: np

Number of function parameters.

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

Function parameters.

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

Generate final summary report.

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

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

logical, intent(in) :: implicit

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

logical, intent(in) :: didvcv

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

logical, intent(in) :: dovcv

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

logical, intent(in) :: redoj

Variable designating whether the Jacobian matrix is to be recomputed for the computation of the covariance matrix (.true.) or not (.false.).

logical, intent(in) :: anajac

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

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

Variable designating why the computations were stopped.

integer, intent(in) :: niter

Number of iterations.

integer, intent(in) :: nfev

Number of function evaluations.

integer, intent(in) :: njev

Number of Jacobian evaluations.

integer, intent(in) :: irank

Rank deficiency of the Jacobian with respect to beta.

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

Approximate reciprocal condition of tfjacb.

integer, intent(in) :: istop

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

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

Sum-of-squares of the weighted epsilons and deltas.

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

Sum-of-squares of the weighted deltas.

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

Sum-of-squares of the weighted epsilons.

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

Penalty parameter for an implicit model.

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

Function parameters.

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

Standard errors of the estimated parameters.

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

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

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

Estimated values of epsilon.

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

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.

public impure subroutine print_errors(info, lunerr, n, m, np, q, ldscld, ldstpd, ldwe, ld2we, ldwd, ld2wd, lrwkmin, liwkmin, fjacb, fjacd, diff, msgb, isodr, msgd, xplusd, nrow, neta, ntol)

Controlling routine for printing error reports.

Arguments

Type IntentOptional Attributes Name
integer, intent(inout) :: info

Variable designating why the computations were stopped.

integer, intent(in) :: lunerr

Logical unit number used for error messages.

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

Minimum acceptable length of array rwork.

integer, intent(in) :: liwkmin

Minimum acceptable length of array iwork.

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

Jacobian with respect to beta.

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

Jacobian with respect to delta.

real(kind=wp), intent(in) :: diff(q,np+m)

Relative differences between the user-supplied and finite difference derivatives for each derivative checked.

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

Error checking results for the Jacobian with respect to beta.

logical, intent(in) :: isodr

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

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

Error checking results for the Jacobian with respect to delta.

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

Values x + delta.

integer, intent(in) :: nrow

Row number of the explanatory variable array at which the derivative is to be checked.

integer, intent(in) :: neta

Number of reliable digits in the model.

integer, intent(in) :: ntol

Number of digits of agreement required between the finite difference and the user-supplied derivatives.

public impure subroutine print_error_inputs(lunerr, info, d1, d2, d3, d4, d5, n, m, q, ldscld, ldstpd, ldwe, ld2we, ldwd, ld2wd, lrwkmin, liwkmin)

Print error reports.

Arguments

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

Logical unit number used for error messages.

integer, intent(inout) :: info

Variable designating why the computations were stopped.

integer, intent(in) :: d1

1st digit (from the left) of info.

integer, intent(in) :: d2

2nd digit (from the left) of info.

integer, intent(in) :: d3

3rd digit (from the left) of info.

integer, intent(in) :: d4

4th digit (from the left) of info.

integer, intent(in) :: d5

5th digit (from the left) of info.

integer, intent(in) :: n

Number of observations.

integer, intent(in) :: m

Number of columns of data in the explanatory variable.

integer, intent(in) :: q

Number of responses per observation.

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

Minimum acceptable length of array rwork.

integer, intent(in) :: liwkmin

Minimum acceptable length of array iwork.

public impure subroutine print_error_derivative(lunerr, n, m, np, q, fjacb, fjacd, diff, msgb1, msgb, isodr, msgd1, msgd, xplusd, nrow, neta, ntol)

Generate the derivative checking report.

Arguments

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

Logical unit number used for error messages.

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.

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

Jacobian with respect to beta.

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

Jacobian with respect to delta.

real(kind=wp), intent(in) :: diff(q,np+m)

Relative differences between the user-supplied and finite difference derivatives for each derivative checked.

integer, intent(in) :: msgb1

Error checking results for the Jacobian with respect to beta.

integer, intent(in) :: msgb(q,np)

Error checking results for the Jacobian with respect to beta.

logical, intent(in) :: isodr

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

integer, intent(in) :: msgd1

Error checking results for the Jacobian with respect to delta.

integer, intent(in) :: msgd(q,m)

Error checking results for the Jacobian with respect to delta.

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

Values of x + delta.

integer, intent(in) :: nrow

Row number of the explanatory variable array at which the derivative is to be checked.

integer, intent(in) :: neta

Number of reliable digits in the model.

integer, intent(in) :: ntol

Number of digits of agreement required between the finite difference and the user-supplied derivatives.

public impure subroutine print_error_fcn(lunerr, d2, d3)

Print error reports indicating that computations were stopped in user-supplied subroutine fcn.

Arguments

Type IntentOptional Attributes Name
integer :: lunerr

Logical unit number used for error messages.

integer :: d2

2nd digit (from the left) of info.

integer :: d3

3rd digit (from the left) of info.