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~doddrv doddrv proc~doddrv->module~odrpack_reports proc~dodmn dodmn proc~dodmn->module~odrpack_reports proc~odr odr proc~odr->module~odrpack_reports

Subroutines

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

Generate initial summary report.

Arguments

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

The value indicating the report to be printed.

integer, intent(in) :: lunrpt

The logical unit number for the computation reports.

logical, intent(in) :: anajac

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

logical, intent(in) :: cdjac

The variable designating whether the Jacobians are computed by central differences (cdjac = .true.) or forward differences (cdjac = .false.).

logical, intent(in) :: chkjac

The variable designating whether the user-supplied Jacobians are to be checked (chkjac = .true.) or not (chkjac = .false.).

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

logical, intent(in) :: restrt

The variable designating whether the call is a restart (restrt = .true.) or not (restrt = .false.).

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

The variable designating whether the covariance matrix is 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.).

integer, intent(in) :: msgb1

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

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

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

integer, intent(in) :: msgd1

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

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

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

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

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

The explanatory variable.

integer, intent(in) :: ldx

The leading dimension of array x.

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

The estimated errors in the explanatory variables.

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.

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

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

The penalty parameter for an implicit model.

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.

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

The scaling values for beta.

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

The lower bounds for beta.

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

The upper bounds for beta.

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. 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

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

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.

public impure subroutine dodpc2(ipr, lunrpt, fstitr, implct, prtpen, pnlty, niter, nfev, wss, actred, prered, alpha, tau, pnorm, np, beta)

Generate iteration reports.

Arguments

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

The value indicating the report to be printed.

integer, intent(in) :: lunrpt

The logical unit number used for computation reports.

logical, intent(in) :: fstitr

The variable designating whether this is the first iteration (fstitr = .true.) or not (fstitr = .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) :: prtpen

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

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

The penalty parameter for an implicit model.

integer, intent(in) :: niter

The number of iterations.

integer, intent(in) :: nfev

The number of function evaluations.

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

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

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

The Levenberg-Marquardt parameter.

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

The trust region diameter.

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

The norm of the scaled estimated parameters.

integer, intent(in) :: np

The number of function parameters.

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

The function parameters.

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)

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.

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)

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.

public impure subroutine dodpe1(lunerr, info, d1, d2, d3, d4, d5, n, m, nq, ldscld, ldstpd, ldwe, ld2we, ldwd, ld2wd, lwkmn, liwkmn)

Print error reports.

Arguments

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

The logical unit number used for error messages.

integer, intent(inout) :: info

The variable designating why the computations were stopped.

integer, intent(in) :: d1

The 1st digit (from the left) of info.

integer, intent(in) :: d2

The 2nd digit (from the left) of info.

integer, intent(in) :: d3

The 3rd digit (from the left) of info.

integer, intent(in) :: d4

The 4th digit (from the left) of info.

integer, intent(in) :: d5

The 5th digit (from the left) of info.

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

The number of responses per observation.

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

The minimum acceptable length of array work.

integer, intent(in) :: liwkmn

The minimum acceptable length of array iwork.

public impure subroutine dodpe2(lunerr, n, m, np, nq, 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

The logical unit number used for error messages.

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.

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

The Jacobian with respect to beta.

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

The Jacobian with respect to delta.

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

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

integer, intent(in) :: msgb1

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

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

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

logical, intent(in) :: isodr

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

integer, intent(in) :: msgd1

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

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

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

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

The values of x + delta.

integer, intent(in) :: nrow

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

integer, intent(in) :: neta

The number of reliable digits in the model.

integer, intent(in) :: ntol

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

public impure subroutine dodpe3(lunerr, d2, d3)

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

Arguments

Type IntentOptional Attributes Name
integer :: lunerr

The logical unit number used for error messages.

integer :: d2

The 2nd digit (from the left) of info.

integer :: d3

The 3rd digit (from the left) of info.

public impure subroutine dodper(info, lunerr, n, m, np, nq, ldscld, ldstpd, ldwe, ld2we, ldwd, ld2wd, lwkmn, liwkmn, 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

The variable designating why the computations were stopped.

integer, intent(in) :: lunerr

The logical unit number used for error messages.

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

The minimum acceptable length of array work.

integer, intent(in) :: liwkmn

The minimum acceptable length of array iwork.

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

The Jacobian with respect to beta.

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

The Jacobian with respect to delta.

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

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

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

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

logical, intent(in) :: isodr

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

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

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

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

The values x + delta.

integer, intent(in) :: nrow

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

integer, intent(in) :: neta

The number of reliable digits in the model.

integer, intent(in) :: ntol

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

public impure subroutine dodphd(head, lunit)

Print report heading.

Arguments

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

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

integer, intent(in) :: lunit

The logical unit number to which the heading is written.