odrpack Module

Main driver routines for finding the weighted explicit or implicit orthogonal distance regression (ODR) or ordinary linear or nonlinear least squares (OLS) solution.


Uses

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

Used by

  • module~~odrpack~~UsedByGraph module~odrpack odrpack module~odrpack_capi odrpack_capi module~odrpack_capi->module~odrpack program~example1 example1 program~example1->module~odrpack program~example2 example2 program~example2->module~odrpack program~example3 example3 program~example3->module~odrpack program~example4 example4 program~example4->module~odrpack program~example5 example5 program~example5->module~odrpack

Subroutines

public impure subroutine odr(fcn, n, m, np, nq, beta, y, x, delta, we, wd, ifixb, ifixx, job, ndigit, taufac, sstol, partol, maxit, iprint, lunerr, lunrpt, stpb, stpd, sclb, scld, work, iwork, info, lower, upper)

Driver routine for finding the weighted explicit or implicit orthogonal distance regression (ODR) or ordinary linear or nonlinear least squares (OLS) solution.

Arguments

Type IntentOptional Attributes Name
procedure(fcn_t) :: fcn

User-supplied subroutine for evaluating the model.

integer, intent(in) :: n

Number of observations.

integer, intent(in) :: m

Number of columns of data in the independent variable.

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: nq

Number of responses per observation.

real(kind=wp), intent(inout) :: beta(:)

Function parameters. Shape: (np).

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

Dependent variable. Shape: (n, nq). Unused when the model is implicit.

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

Explanatory variable. Shape: (n, m).

real(kind=wp), intent(inout), optional :: delta(:,:)

Error in the x data. Shape: (n, m). Initial guess on input and estimated value on output.

real(kind=wp), intent(in), optional :: we(:,:,:)

epsilon weights. Shape: (1<=ldwe<=n, 1<=ld2we<=nq, nq). See p. 25.

real(kind=wp), intent(in), optional :: wd(:,:,:)

delta weights. Shape: (1<=ldwd<=n, 1<=ld2wd<=m, m). See p. 26.

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

Values designating whether the elements of beta are fixed at their input values or not. Shape: (np).

integer, intent(in), optional :: ifixx(:,:)

Values designating whether the elements of x are fixed at their input values or not. Shape: (1<=ldifx<=n, m). See p. 27.

integer, intent(in), optional :: job

Variable controlling problem initialization and computational method.

integer, intent(in), optional :: ndigit

Number of accurate digits in the function results, as supplied by the user.

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

Factor used to compute the initial trust region diameter.

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

Sum-of-squares convergence stopping tolerance.

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

Parameter convergence stopping tolerance.

integer, intent(in), optional :: maxit

Maximum number of iterations allowed.

integer, intent(in), optional :: iprint

Print control variable.

integer, intent(in), optional :: lunerr

Logical unit number for error messages. Available options are: 0 => no output. 6 => output to standard error. other => output to logical unit number lunerr.

integer, intent(in), optional :: lunrpt

Logical unit number for computation reports. Available options are: 0 => no output. 6 => output to standard error. other => output to logical unit number lunrpt.

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

Relative step for computing finite difference derivatives with respect to beta. Shape: (np).

real(kind=wp), intent(in), optional :: stpd(:,:)

Relative step for computing finite difference derivatives with respect to delta. Shape: (1<=ldstpd<=n, m). See p. 31.

real(kind=wp), intent(in), optional :: sclb(:)

Scaling values for beta. Shape: (np).

real(kind=wp), intent(in), optional :: scld(:,:)

Scaling values for delta. Shape: (1<=ldscld<=n, m). See p. 32.

real(kind=wp), intent(inout), optional, target :: work(:)

Real work space.

integer, intent(inout), optional, target :: iwork(:)

Integer work space.

integer, intent(out), optional :: info

Variable designating why the computations were stopped.

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

Lower bound on beta. Shape: (np).

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

Upper bound on beta. Shape: (np).

public impure subroutine dodcnt(fcn, n, m, np, nq, beta, y, ldy, x, ldx, we, ldwe, ld2we, wd, ldwd, ld2wd, ifixb, ifixx, ldifx, job, ndigit, taufac, sstol, partol, maxit, iprint, lunerr, lunrpt, stpb, stpd, ldstpd, sclb, scld, ldscld, work, lwork, tempret, iwork, liwork, info, lower, upper)

Driver routine for finding the weighted explicit or implicit orthogonal distance regression (ODR) or ordinary linear or nonlinear least squares (OLS) solution.

Arguments

Type IntentOptional Attributes Name
procedure(fcn_t) :: fcn

User-supplied subroutine for evaluating the model.

integer, intent(in) :: n

The number of observations.

integer, intent(in) :: m

The number of columns of data in the independent variable.

integer, intent(in) :: np

The number of function parameters.

integer, intent(in) :: nq

The number of responses per observation.

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

The function parameters.

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

The dependent 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 independent variable.

integer, intent(in) :: ldx

The leading dimension of array x.

real(kind=wp), intent(inout) :: 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.

integer, intent(inout) :: job

The variable controlling problem initialization and computational method.

integer, intent(in) :: ndigit

The number of accurate digits in the function results, as 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 user-supplied parameter convergence stopping tolerance.

integer, intent(in) :: maxit

The maximum number of iterations allowed.

integer, intent(in) :: iprint

The print control variables.

integer, intent(in) :: lunerr

The logical unit number used for error messages.

integer, intent(in) :: lunrpt

The logical unit number used for computation reports.

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

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

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

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

integer, intent(in) :: ldstpd

The leading dimension of array stpd.

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

The scaling values for beta.

real(kind=wp), intent(in) :: scld(ldscld,m)

The scaling value for delta.

integer, intent(in) :: ldscld

The leading dimension of array scld.

real(kind=wp), intent(inout) :: work(lwork)

The real work space.

integer, intent(in) :: lwork

The length of vector work.

real(kind=wp), intent(inout) :: tempret(:,:)

Temporary work array for holding return values before copying to a lower rank array.

integer, intent(inout) :: iwork(liwork)

The integer work space.

integer, intent(in) :: liwork

The length of vector iwork.

integer, intent(out) :: info

The variable designating why the computations were stopped.

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

The lower bound on beta.

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

The upper bound on beta.

public impure subroutine doddrv(head, fstitr, prtpen, fcn, n, m, np, nq, beta, y, ldy, x, ldx, we, ldwe, ld2we, wd, ldwd, ld2wd, ifixb, ifixx, ldifx, job, ndigit, taufac, sstol, partol, maxit, iprint, lunerr, lunrpt, stpb, stpd, ldstpd, sclb, scld, ldscld, work, lwork, tempret, iwork, liwork, maxit1, tstimp, info, lower, upper)

Perform error checking and initialization, and begin procedure for performing orthogonal distance regression (ODR) or ordinary linear or nonlinear least squares (OLS).

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

logical, intent(inout) :: fstitr

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

logical, intent(inout) :: prtpen

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

procedure(fcn_t) :: fcn

The user-supplied subroutine for evaluating the model.

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

The function parameters.

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

The dependent 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(inout) :: 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.

integer, intent(inout) :: job

The variable controlling problem initialization and computational method.

integer, intent(in) :: ndigit

The number of accurate digits in the function results, as 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.

integer, intent(in) :: iprint

The print control variable.

integer, intent(in) :: lunerr

The logical unit number used for error messages.

integer, intent(in) :: lunrpt

The logical unit number used for computation reports.

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

The step size for finite difference derivatives with respect to beta.

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

The step size for finite difference derivatives with respect to delta.

integer, intent(in) :: ldstpd

The leading dimension of array stpd.

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

The scaling values for beta.

real(kind=wp), intent(in) :: scld(ldscld,m)

The scaling values for delta.

integer, intent(in) :: ldscld

The leading dimension of array scld.

real(kind=wp), intent(inout) :: work(lwork)

The real work space.

integer, intent(in) :: lwork

The length of vector work.

real(kind=wp), intent(inout) :: tempret(:,:)

Temporary work array for holding return values before copying to a lower rank array.

integer, intent(inout) :: iwork(liwork)

The integer work space.

integer, intent(in) :: liwork

The length of vector iwork.

integer, intent(out) :: maxit1

For implicit models, the iterations allowed for the next penalty parameter value.

real(kind=wp), intent(out) :: tstimp

The relative change in the parameters between the initial values and the solution.

integer, intent(inout) :: info

The variable designating why the computations were stopped.

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

The lower bound for beta.

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

The upper bound for beta.

public impure subroutine dodmn(head, fstitr, prtpen, fcn, n, m, np, nq, job, beta, y, ldy, x, ldx, we, we1, ldwe, ld2we, wd, ldwd, ld2wd, ifixb, ifixx, ldifx, betac, betan, betas, s, delta, deltan, deltas, lower, upper, t, f, fn, fs, fjacb, msgb, fjacd, msgd, ssf, ss, tt, ldtt, stpb, stpd, ldstpd, xplusd, wrk, lwrk, work, lwork, tempret, iwork, liwork, info, bound)

Iteratively compute least squares solution.

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

logical, intent(inout) :: fstitr

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

logical, intent(inout) :: prtpen

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

procedure(fcn_t) :: fcn

The user supplied subroutine for evaluating the model.

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(inout) :: job

The variable controlling problem initialization and computational method.

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

The function parameters.

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

The dependent 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) :: we(ldwe,ld2we,nq)

The epsilon weights.

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

The square root of the epsilon weights.

integer, intent(in) :: ldwe

The leading dimension of arrays we and we1.

integer, intent(in) :: ld2we

The second dimension of arrays we and we1.

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

The current estimated values of the unfixed betas.

real(kind=wp), intent(out) :: betan(np)

The new estimated values of the unfixed betas.

real(kind=wp), intent(inout) :: betas(np)

The saved estimated values of the unfixed betas.

real(kind=wp), intent(out) :: s(np)

The step for beta.

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

The estimated errors in the explanatory variables.

real(kind=wp), intent(out) :: deltan(n,m)

The new estimated errors in the explanatory variables.

real(kind=wp), intent(inout) :: deltas(n,m)

The saved estimated errors in the explanatory variables.

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

The lower bound for unfixed betas.

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

The upper bound for unfixed betas.

real(kind=wp), intent(out) :: t(n,m)

The step for delta.

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

The (weighted) estimated values of epsilon.

real(kind=wp), intent(out) :: fn(n,nq)

The new predicted values from the function.

real(kind=wp), intent(out) :: fs(n,nq)

The saved predicted values from the function.

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

The Jacobian with respect to beta.

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

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

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

The Jacobian with respect to delta.

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

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

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

The scaling values used for beta.

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

The scaling values used for the unfixed betas.

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

The scaling values used 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 each 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.

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

The values of x + delta.

real(kind=wp), intent(inout) :: wrk(lwrk)

A work array, equivalenced to wrk1 and wrk2.

integer, intent(in) :: lwrk

The length of vector wrk.

real(kind=wp), intent(inout) :: work(lwork)

The real (wp) workspace.

integer, intent(in) :: lwork

The length of vector work.

real(kind=wp), intent(inout) :: tempret(:,:)

Temporary work array for holding return values before copying to a lower rank array.

integer, intent(inout) :: iwork(liwork)

The integer workspace.

integer, intent(in) :: liwork

The length of vector iwork.

integer, intent(inout) :: info

The variable designating why the computations were stopped.

integer, intent(out) :: bound(np)

The values of the bounds for beta.

public pure subroutine workspace_dimensions(n, m, np, nq, isodr, lwork, liwork)

Calculate the dimensions of the workspace arrays.

Arguments

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

Number of observations.

integer, intent(in) :: m

Number of columns of data in the independent variable.

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: nq

Number of responses per observation.

logical, intent(in) :: isodr

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

integer, intent(out) :: lwork

Length of work array.

integer, intent(out) :: liwork

Length of iwork array.