odrpack_core Module

Core mathematical routines, except drivers, and BLAS/LINPACK.


Uses

  • module~~odrpack_core~~UsesGraph module~odrpack_core 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

Used by

  • module~~odrpack_core~~UsedByGraph module~odrpack_core odrpack_core module~odrpack_capi odrpack_capi module~odrpack_capi->module~odrpack_core proc~dodcnt dodcnt proc~dodcnt->module~odrpack_core proc~doddrv doddrv proc~doddrv->module~odrpack_core proc~dodmn dodmn proc~dodmn->module~odrpack_core proc~dodpc1 dodpc1 proc~dodpc1->module~odrpack_core proc~dodpc3 dodpc3 proc~dodpc3->module~odrpack_core proc~dodpcr dodpcr proc~dodpcr->module~odrpack_core proc~odr odr proc~odr->module~odrpack_core

Abstract Interfaces

abstract interface

  • public subroutine fcn_t(n, m, np, nq, ldn, ldm, ldnp, beta, xplusd, ifixb, ifixx, ldifx, ideval, f, fjacb, fjacd, istop)

    User-supplied subroutine for evaluating the model.

    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.

    integer, intent(in) :: ldn

    Leading dimension declarator equal or exceeding n.

    integer, intent(in) :: ldm

    Leading dimension declarator equal or exceeding m.

    integer, intent(in) :: ldnp

    Leading dimension declarator equal or exceeding np.

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

    Current values of parameters.

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

    Current value of explanatory variable, i.e., x + delta.

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

    Indicators for "fixing" parameters (beta).

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

    Indicators for "fixing" explanatory variable (x).

    integer, intent(in) :: ldifx

    Leading dimension of array ifixx.

    integer, intent(in) :: ideval

    Indicator for selecting computation to be performed.

    real(kind=wp), intent(out) :: f(ldn,nq)

    Predicted function values.

    real(kind=wp), intent(out) :: fjacb(ldn,ldnp,nq)

    Jacobian with respect to beta.

    real(kind=wp), intent(out) :: fjacd(ldn,ldm,nq)

    Jacobian with respect to errors delta.

    integer, intent(out) :: istop

    Stopping condition, with meaning as follows. 0 means current beta and x+delta were acceptable and values were computed successfully. 1 means current beta and x+delta are not acceptable; 'odrpack' should select values closer to most recently used values if possible. -1 means current beta and x+delta are not acceptable; 'odrpack' should stop.


Functions

public pure function derstep(itype, k, betak, ssf, stpb, neta) result(derstepr)

Compute step size for center and forward difference calculations.

Arguments

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

The finite difference method being used, where: itype = 0 indicates forward finite differences, and itype = 1 indicates central finite differences.

integer, intent(in) :: k

Index into beta where betak resides.

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

The k-th function parameter.

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

The scale used for the betas.

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

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

integer, intent(in) :: neta

Number of good digits in the function results.

Return Value real(kind=wp)

public pure function dhstep(itype, neta, i, j, stp, ldstp) result(dhstepr)

Set relative step size for finite difference derivatives.

Arguments

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

The finite difference method being used, where: itype = 0 indicates forward finite differences, and itype = 1 indicates central finite differences.

integer, intent(in) :: neta

The number of good digits in the function results.

integer, intent(in) :: i

An identifier for selecting user-supplied step sizes.

integer, intent(in) :: j

An identifier for selecting user-supplied step sizes.

real(kind=wp), intent(in) :: stp(ldstp,j)

The step size for the finite difference derivative.

integer, intent(in) :: ldstp

The leading dimension of array stp.

Return Value real(kind=wp)

public pure function dppnml(p) result(dppnmlr)

Compute the percent point function value for the normal (Gaussian) distribution with mean 0 and standard deviation 1, and with probability density function:

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: p

The probability at which the percent point is to be evaluated. p must lie between 0.0 and 1.0, exclusive.

Return Value real(kind=wp)

public pure function dppt(p, idf) result(dpptr)

Compute the percent point function value for the student's T distribution with idf degrees of freedom.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: p

The probability at which the percent point is to be evaluated. p must lie between 0.0 and 1.0, exclusive.

integer, intent(in) :: idf

The (positive integer) degrees of freedom.

Return Value real(kind=wp)


Subroutines

public subroutine dodlm(n, m, np, nq, npp, f, fjacb, fjacd, wd, ldwd, ld2wd, ss, tt, ldtt, delta, alpha2, tau, epsfcn, isodr, tfjacb, omega, u, qraux, jpvt, s, t, nlms, rcond, irank, wrk1, wrk2, wrk3, wrk4, wrk5, wrk, lwrk, tempret, istopc)

Compute Levenberg-Marquardt parameter and steps s and t using analog of the trust-region Levenberg-Marquardt algorithm.

Arguments

Type IntentOptional Attributes Name
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.

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

The (weighted) estimated values of epsilon.

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

The scaling values used for the unfixed betas.

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

The scale used for the deltas.

integer, intent(in) :: ldtt

The leading dimension of array tt.

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

The estimated errors in the explanatory variables.

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

The current Levenberg-Marquardt parameter.

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

The trust region diameter.

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

The function's precision.

logical, intent(in) :: isodr

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

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

The array omega*fjacb.

real(kind=wp), intent(out) :: omega(nq,nq)

The array (I-fjacd*inv(p)*trans(fjacd))**(-1/2).

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

The approximate null vector for tfjacb.

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

The array required to recover the orthogonal part of the Q-R decomposition.

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

The pivot vector.

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

The step for beta.

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

The step for delta.

integer, intent(out) :: nlms

The number of Levenberg-Marquardt steps taken.

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

The approximate reciprocal condition of tfjacb.

integer, intent(out) :: irank

The rank deficiency of the Jacobian wrt beta.

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

A work array of (n, nq, m) elements.

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

A work array of (n, nq) elements.

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

A work array of (np) elements.

real(kind=wp), intent(out) :: wrk4(m,m)

A work array of (m, m) elements.

real(kind=wp), intent(out) :: wrk5(m)

A work array of (m) elements.

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

A work array of (lwrk) elements, equivalenced to wrk1 and wrk2.

integer, intent(in) :: lwrk

The length of vector wrk.

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

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

integer, intent(out) :: istopc

The variable designating whether the computations were stopped due to some other numerical error detected within subroutine dodstp.

public pure subroutine dacces(n, m, np, nq, ldwe, ld2we, work, lwork, iwork, liwork, access, isodr, jpvt, omega, u, qraux, sd, vcv, wrk1, wrk2, wrk3, wrk4, wrk5, wrk6, nnzw, npp, job, partol, sstol, maxit, taufac, eta, neta, lunrpt, ipr1, ipr2, ipr2f, ipr3, wss, rvar, idf, tau, alpha, niter, nfev, njev, int2, olmavg, rcond, irank, actrs, pnorm, prers, rnorms, istop)

Access or store values in the work arrays.

Arguments

Type IntentOptional Attributes Name
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) :: ldwe

The leading dimension of array we.

integer, intent(in) :: ld2we

The second dimension of array we.

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

The real work space.

integer, intent(in) :: lwork

The length of vector work.

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

The integer work space.

integer, intent(in) :: liwork

The length of vector iwork.

logical, intent(in) :: access

The variable designating whether information is to be accessed from the work arrays (access = .true.) or stored in them (access = .false.).

logical, intent(in) :: isodr

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

integer, intent(out) :: jpvt

The pivot vector.

integer, intent(out) :: omega

The starting location in array work of array omega.

integer, intent(out) :: u

The starting location in array work of array u.

integer, intent(out) :: qraux

The starting location in array work of array qraux.

integer, intent(out) :: sd

The starting location in array work of array sd.

integer, intent(out) :: vcv

The starting location in array work of array vcv.

integer, intent(out) :: wrk1

The starting location in array work of array wrk1.

integer, intent(out) :: wrk2

The starting location in array work of array wrk2.

integer, intent(out) :: wrk3

The starting location in array work of array wrk3.

integer, intent(out) :: wrk4

The starting location in array work of array wrk4.

integer, intent(out) :: wrk5

The starting location in array work of array wrk5.

integer, intent(out) :: wrk6

The starting location in array work of array wrk6.

integer, intent(out) :: nnzw

The number of nonzero weighted observations.

integer, intent(out) :: npp

The number of function parameters actually estimated.

integer, intent(out) :: job

The variable controlling problem initialization and computational method.

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

The parameter convergence stopping tolerance.

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

The sum-of-squares convergence stopping tolerance.

integer, intent(out) :: maxit

The maximum number of iterations allowed.

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

The factor used to compute the initial trust region diameter.

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

The relative noise in the function results.

integer, intent(out) :: neta

The number of accurate digits in the function results.

integer, intent(out) :: lunrpt

The logical unit number used for computation reports.

integer, intent(out) :: ipr1

The value of the fourth digit (from the right) of iprint, which controls the initial summary report.

integer, intent(out) :: ipr2

The value of the third digit (from the right) of iprint, which controls the iteration reports.

integer, intent(out) :: ipr2f

The value of the second digit (from the right) of iprint, which controls the frequency of the iteration reports.

integer, intent(out) :: ipr3

The value of the first digit (from the right) of iprint, which controls the final summary report.

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

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

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

The residual variance, i.e. the standard deviation squared.

integer, intent(inout) :: 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(inout) :: tau

The trust region diameter.

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

The Levenberg-Marquardt parameter.

integer, intent(inout) :: niter

The number of iterations taken.

integer, intent(inout) :: nfev

The number of function evaluations.

integer, intent(inout) :: njev

The number of Jacobian evaluations.

integer, intent(inout) :: int2

The number of internal doubling steps.

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

The average number of Levenberg-Marquardt steps per iteration.

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

The approximate reciprocal condition of fjacb.

integer, intent(inout) :: irank

The rank deficiency of the Jacobian wrt beta.

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

The saved actual relative reduction in the sum-of-squares.

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

The norm of the scaled estimated parameters.

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

The saved predicted relative reduction in the sum-of-squares.

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

The norm of the saved weighted epsilons and deltas.

integer, intent(inout) :: istop

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

public pure subroutine desubi(n, m, wd, ldwd, ld2wd, alpha, tt, ldtt, i, e)

Compute e = wd + alpha*tt**2.

Arguments

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

The number of observations.

integer, intent(in) :: m

The number of columns of data in the independent variable.

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

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

The Levenberg-Marquardt parameter.

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

The scaling values used for delta.

integer, intent(in) :: ldtt

The leading dimension of array tt.

integer, intent(in) :: i

An indexing variable.

real(kind=wp), intent(out) :: e(m,m)

The value of the array e = wd + alpha*tt**2.

public subroutine detaf(fcn, n, m, np, nq, xplusd, beta, epsmac, nrow, partmp, pv0, ifixb, ifixx, ldifx, istop, nfev, eta, neta, wrk1, wrk2, wrk6, wrk7, info, lower, upper)

Compute noise and number of good digits in function results.

Arguments

Type IntentOptional Attributes Name
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(in) :: xplusd(n,m)

The values of x + delta.

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

The function parameters.

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

The value of machine precision.

integer, intent(in) :: nrow

The row number at which the derivative is to be checked.

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

The model parameters.

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

The original predicted values.

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(out) :: istop

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

integer, intent(inout) :: nfev

The number of function evaluations.

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

The noise in the model results.

integer, intent(out) :: neta

The number of accurate digits in the model results.

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

A work array of (n, m, nq) elements.

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

A work array of (n, nq) elements.

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

A work array of (n, np, nq) elements.

real(kind=wp), intent(out) :: wrk7(-2:2,nq)

A work array of (5, nq) elements.

integer, intent(out) :: info

The variable indicating the status of the computation.

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

The lower bound of beta.

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

The upper bound of beta.

public subroutine devjac(fcn, anajac, cdjac, n, m, np, nq, betac, beta, stpb, ifixb, ifixx, ldifx, x, ldx, delta, xplusd, stpd, ldstpd, ssf, tt, ldtt, neta, fn, stp, wrk1, wrk2, wrk3, wrk6, tempret, fjacb, isodr, fjacd, we1, ldwe, ld2we, njev, nfev, istop, info, lower, upper)

Compute the weighted Jacobians wrt beta and delta.

Arguments

Type IntentOptional Attributes Name
procedure(fcn_t) :: fcn

The user-supplied subroutine for evaluating the model.

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 by forward differences (cdjac = .false.).

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

The current estimated values of the unfixed betas.

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

The function parameters.

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

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

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 delta are fixed at their input values or not.

integer, intent(in) :: ldifx

The leading dimension of array ifixx.

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

The estimated values of delta.

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

The values of x + delta.

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

The scale used for the 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.

integer, intent(in) :: neta

The number of accurate digits in the function results.

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

The predicted values of the function at the current point.

real(kind=wp), intent(out) :: stp(n)

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

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

A work array of (n, m, nq) elements.

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

A work array of (n, nq) elements.

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

A work array of (np) elements.

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

A work array of (n, np, nq) elements.

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

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

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

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

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

The Jacobian with respect to delta.

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

The square roots of the epsilon weights in array we.

integer, intent(in) :: ldwe

The leading dimension of arrays we and we1.

integer, intent(in) :: ld2we

The second dimension of arrays we and we1.

integer, intent(inout) :: njev

The number of Jacobian evaluations.

integer, intent(inout) :: nfev

The number of function evaluations.

integer, intent(out) :: istop

The variable designating that the user wishes the computations stopped.

integer, intent(out) :: info

The variable designating why the computations were stopped.

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

The lower bound of beta.

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

The upper bound of beta.

public pure subroutine dfctr(oksemi, a, lda, n, info)

Factor the positive (semi)definite matrix a using a modified Cholesky factorization.

Arguments

Type IntentOptional Attributes Name
logical, intent(in) :: oksemi

The indicating whether the factored array can be positive semidefinite (oksemi = .true.) or whether it must be found to be positive definite (oksemi = .false.).

real(kind=wp), intent(inout) :: a(lda,n)

The array to be factored. Upon return, a contains the upper triangular matrix r so that a = trans(r)*r where the strict lower triangle is set to zero. If info /= 0, the factorization is not complete.

integer, intent(in) :: lda

The leading dimension of array a.

integer, intent(in) :: n

The number of rows and columns of data in array a.

integer, intent(out) :: info

An indicator variable, where if: info = 0 then factorization was completed. info = k signals an error condition. The leading minor of order k is not positive (semi)definite.

public pure subroutine dfctrw(n, m, nq, npp, isodr, we, ldwe, ld2we, wd, ldwd, ld2wd, wrk0, wrk4, we1, nnzw, info)

Check input parameters, indicating errors found using nonzero values of argument info as described in the ODRPACK95 reference guide.

Arguments

Type IntentOptional Attributes Name
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) :: npp

The number of function parameters being estimated.

logical, intent(in) :: isodr

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

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

The (squared) 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 (squared) delta weights.

integer, intent(in) :: ldwd
integer, intent(in) :: ld2wd

The second dimension of array wd.

real(kind=wp), intent(out) :: wrk0(nq,nq)

A work array of (nq, nq) elements.

real(kind=wp), intent(out) :: wrk4(m,m)

A work array of (m, m) elements.

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

The factored epsilon weights, such that trans(we1)*we1 = we.

integer, intent(out) :: nnzw

The number of nonzero weighted observations.

integer, intent(out) :: info

The variable designating why the computations were stopped.

public pure subroutine dflags(job, restrt, initd, dovcv, redoj, anajac, cdjac, chkjac, isodr, implct)

Set flags indicating conditions specified by job.

Arguments

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

The variable controlling problem initialization and computational method.

logical, intent(out) :: restrt

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

logical, intent(out) :: initd

The variable designating whether delta is to be initialized to zero (initd = .true.) or to the first n by m elements of array work (initd = .false.).

logical, intent(out) :: dovcv

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

logical, intent(out) :: 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(out) :: anajac

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

logical, intent(out) :: cdjac

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

logical, intent(out) :: chkjac

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

logical, intent(out) :: isodr

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

logical, intent(out) :: implct

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

public pure subroutine difix(n, m, ifix, ldifix, t, ldt, tfix, ldtfix)

Set elements of t to zero according to ifix.

Arguments

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

The number of rows of data in the array.

integer, intent(in) :: m

The number of columns of data in the array.

integer, intent(in) :: ifix(ldifix,m)

The array designating whether an element of t is to be set to zero.

integer, intent(in) :: ldifix

The leading dimension of array ifix.

real(kind=wp), intent(in) :: t(ldt,m)

The array being set to zero according to the elements of ifix.

integer, intent(in) :: ldt

The leading dimension of array t.

real(kind=wp), intent(out) :: tfix(ldtfix,m)

The resulting array.

integer, intent(in) :: ldtfix

The leading dimension of array tfix.

public pure subroutine diwinf(m, np, nq, msgbi, msgdi, ifix2i, istopi, nnzwi, nppi, idfi, jobi, iprini, luneri, lunrpi, nrowi, ntoli, netai, maxiti, niteri, nfevi, njevi, int2i, iranki, ldtti, boundi, liwkmn)

Get storage locations within integer work space.

Arguments

Type IntentOptional Attributes Name
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.

integer, intent(out) :: msgbi

The starting location in array iwork of array msgb.

integer, intent(out) :: msgdi

The starting location in array iwork of array msgd.

integer, intent(out) :: ifix2i

The starting location in array iwork of array ifix2.

integer, intent(out) :: istopi

The location in array iwork of variable istop.

integer, intent(out) :: nnzwi

The location in array iwork of variable nnzw.

integer, intent(out) :: nppi

The location in array iwork of variable npp.

integer, intent(out) :: idfi

The location in array iwork of variable idf.

integer, intent(out) :: jobi

The location in array iwork of variable job.

integer, intent(out) :: iprini

The location in array iwork of variable iprint.

integer, intent(out) :: luneri

The location in array iwork of variable lunerr.

integer, intent(out) :: lunrpi

The location in array iwork of variable lunrpt.

integer, intent(out) :: nrowi

The location in array iwork of variable nrow.

integer, intent(out) :: ntoli

The location in array iwork of variable ntol.

integer, intent(out) :: netai

The location in array iwork of variable neta.

integer, intent(out) :: maxiti

The location in array iwork of variable maxit.

integer, intent(out) :: niteri

The location in array iwork of variable niter.

integer, intent(out) :: nfevi

The location in array iwork of variable nfev.

integer, intent(out) :: njevi

The location in array iwork of variable njev.

integer, intent(out) :: int2i

The location in array iwork of variable int2.

integer, intent(out) :: iranki

The location in array iwork of variable irank.

integer, intent(out) :: ldtti

The location in array iwork of variable ldtt.

integer, intent(out) :: boundi

The location in array iwork of variable bound.

integer, intent(out) :: liwkmn

The minimum acceptable length of array iwork.

public pure subroutine diniwk(n, m, np, work, lwork, iwork, liwork, x, ldx, ifixx, ldifx, scld, ldscld, beta, sclb, sstol, partol, maxit, taufac, job, iprint, lunerr, lunrpt, lower, upper, epsmai, sstoli, partli, maxiti, taufci, jobi, iprini, luneri, lunrpi, ssfi, tti, ldtti, deltai, loweri, upperi, boundi)

Initialize work vectors as necessary.

Arguments

Type IntentOptional Attributes Name
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.

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

The real work space.

integer, intent(in) :: lwork

The length of vector work.

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

The integer work space.

integer, intent(in) :: liwork

The length of vector iwork.

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

The independent 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) :: scld(ldscld,m)

The scaling values for delta.

integer, intent(in) :: ldscld

The leading dimension of array scld.

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

The function parameters.

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

The scaling values for beta.

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

The sum-of-squares convergence stopping criteria.

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

The parameter convergence stopping criteria.

integer, intent(in) :: maxit

The maximum number of iterations allowed.

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

The factor used to compute the initial trust region diameter.

integer, intent(in) :: job

The variable controlling problem initialization and computational method.

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

The lower bounds for the function parameters.

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

The upper bounds for the function parameters.

integer, intent(in) :: epsmai

The location in array work of variable epsmac.

integer, intent(in) :: sstoli

The location in array work of variable sstol.

integer, intent(in) :: partli

The location in array work of variable partol.

integer, intent(in) :: maxiti

The location in array iwork of variable maxit.

integer, intent(in) :: taufci

The location in array work of variable taufac.

integer, intent(in) :: jobi

The location in array iwork of variable job.

integer, intent(in) :: iprini

The location in array iwork of variable iprint.

integer, intent(in) :: luneri

The location in array iwork of variable lunerr.

integer, intent(in) :: lunrpi

The location in array iwork of variable lunrpt.

integer, intent(in) :: ssfi

The starting location in array work of array ssf.

integer, intent(in) :: tti

The starting location in array work of the array tt.

integer, intent(in) :: ldtti

The leading dimension of array tt.

integer, intent(in) :: deltai

The starting location in array work of array delta.

integer, intent(in) :: loweri

The starting location in array iwork of array lower.

integer, intent(in) :: upperi

The starting location in array iwork of array upper.

integer, intent(in) :: boundi

The location in array iwork of variable bound.

public subroutine djaccd(fcn, n, m, np, nq, beta, x, ldx, delta, xplusd, ifixb, ifixx, ldifx, stpb, stpd, ldstpd, ssf, tt, ldtt, neta, fn, stp, wrk1, wrk2, wrk3, wrk6, fjacb, isodr, fjacd, nfev, istop, info, lower, upper)

Compute central difference approximations to the Jacobian wrt the estimated betas and wrt the deltas.

Arguments

Type IntentOptional Attributes Name
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) :: 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(inout) :: xplusd(n,m)

The values of x + delta.

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

integer, intent(in) :: ldstpd

The leading dimension of array stpd.

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

The scaling values used for beta.

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

The scaling values used for delta.

integer, intent(in) :: ldtt

The leading dimension of array tt.

integer, intent(in) :: neta

The number of good digits in the function results.

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

The new predicted values from the function. Used when parameter is on a boundary.

real(kind=wp), intent(out) :: stp(n)

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

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

A work array of (n, m, nq) elements.

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

A work array of (n, nq) elements.

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

A work array of (np) elements.

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

A work array of (n, np, nq) elements.

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

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

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

The Jacobian with respect to delta.

integer, intent(inout) :: nfev

The number of function evaluations.

integer, intent(out) :: istop

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

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 subroutine djacfd(fcn, n, m, np, nq, beta, x, ldx, delta, xplusd, ifixb, ifixx, ldifx, stpb, stpd, ldstpd, ssf, tt, ldtt, neta, fn, stp, wrk1, wrk2, wrk3, wrk6, fjacb, isodr, fjacd, nfev, istop, info, lower, upper)

Compute forward difference approximations to the Jacobian wrt the estimated betas and wrt the deltas.

Arguments

Type IntentOptional Attributes Name
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) :: 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(inout) :: xplusd(n,m)

The values of x + delta.

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

integer, intent(in) :: ldstpd

The leading dimension of array stpd.

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

The scaling values used for beta.

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

The scaling values used for delta.

integer, intent(in) :: ldtt

The leading dimension of array tt.

integer, intent(in) :: neta

The number of good digits in the function results.

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

The new predicted values from the function. Used when parameter is on a boundary.

real(kind=wp), intent(out) :: stp(n)

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

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

A work array of (n, m, nq) elements.

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

A work array of (n, nq) elements.

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

A work array of (np) elements.

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

A work array of (n, np, nq) elements.

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

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

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

The Jacobian with respect to delta.

integer, intent(inout) :: nfev

The number of function evaluations.

integer, intent(out) :: istop

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

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 subroutine djck(fcn, n, m, np, nq, beta, betaj, xplusd, ifixb, ifixx, ldifx, stpb, stpd, ldstpd, ssf, tt, ldtt, eta, neta, ntol, nrow, isodr, epsmac, pv0i, fjacb, fjacd, msgb, msgd, diff, istop, nfev, njev, wrk1, wrk2, wrk6, interval)

Driver routine for the derivative checking process.

Arguments

Type IntentOptional Attributes Name
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(inout) :: betaj(np)

The function parameters offset such that steps don't cross bounds.

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

The values of x + delta.

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

The step size for finite difference derivatives wrt beta.

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

The step size for finite difference derivatives wrt delta.

integer, intent(in) :: ldstpd

The leading dimension of array stpd.

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

The scaling values used for beta.

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

The relative noise in the function results.

integer, intent(in) :: neta

The number of reliable digits in the model results.

integer, intent(out) :: ntol

The number of digits of agreement required between the numerical derivatives and the user supplied derivatives.

integer, intent(in) :: nrow

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

logical, intent(in) :: isodr

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

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

The value of machine precision.

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

The predicted values using the user supplied parameter estimates.

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

The Jacobian with respect to beta.

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

The Jacobian with respect to delta.

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

The error checking results for the Jacobian wrt beta.

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

The error checking results for the Jacobian wrt delta.

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

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

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

The number of function evaluations.

integer, intent(inout) :: njev

The number of Jacobian evaluations.

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

A work array of (n, m, nq) elements.

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

A work array of (n, nq) elements.

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

A work array of (n, np, nq) elements.

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

Specifies which checks can be performed when checking derivatives based on the interval of the bound constraints.

public subroutine djckc(fcn, n, m, np, nq, beta, xplusd, ifixb, ifixx, ldifx, eta, tol, nrow, epsmac, j, lq, hc, iswrtb, fd, typj, pvpstp, stp0, pv, d, diffj, msg, istop, nfev, wrk1, wrk2, wrk6)

Check whether high curvature could be the cause of the disagreement between the numerical and analytic derviatives.

Arguments

Type IntentOptional Attributes Name
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(inout) :: xplusd(n,m)

The values of x + delta.

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

The relative noise in the model.

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

The agreement tolerance.

integer, intent(in) :: nrow

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

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

The value of machine precision.

integer, intent(in) :: j

The index of the partial derivative being examined.

integer, intent(in) :: lq

The response currently being examined.

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

The relative step size for central finite differences.

logical, intent(in) :: iswrtb

The variable designating whether the derivatives wrt beta (iswrtb = .true.) or delta (iswrtb = .false.) are being checked.

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

The forward difference derivative wrt the j-th parameter.

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

The typical size of the j-th unknown beta or delta.

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

The predicted value for row nrow of the model based on the current parameter estimates for all but the j-th parameter value, which is beta(j) + stp0.

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

The initial step size for the finite difference derivative.

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

The predicted value of the model for row nrow.

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

The derivative with respect to the j-th unknown parameter.

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

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

integer, intent(out) :: msg(nq,j)

The error checking results.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

The number of function evaluations.

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

A work array of (n, m, nq) elements.

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

A work array of (n, nq) elements.

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

A work array of (n, np, nq) elements.

public subroutine djckf(fcn, n, m, np, nq, beta, xplusd, ifixb, ifixx, ldifx, eta, tol, nrow, j, lq, iswrtb, fd, typj, pvpstp, stp0, curve, pv, d, diffj, msg, istop, nfev, wrk1, wrk2, wrk6)

Check whether finite precision arithmetic could be the cause of the disagreement between the derivatives.

Arguments

Type IntentOptional Attributes Name
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(inout) :: xplusd(n,m)

The values of x + delta.

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

The relative noise in the model.

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

The agreement tolerance.

integer, intent(in) :: nrow

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

integer, intent(in) :: j

The index of the partial derivative being examined.

integer, intent(in) :: lq

The response currently being examined.

logical, intent(in) :: iswrtb

The variable designating whether the derivatives wrt beta (iswrtb = .true.) or delta (iswrtb = .false.) are being checked.

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

The forward difference derivative wrt the j-th parameter.

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

The typical size of the j-th unknown beta or delta.

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

The predicted value for row nrow of the model based on the current parameter estimates for all but the j-th parameter value, which is beta(j) + stp0.

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

The step size for the finite difference derivative.

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

A measure of the curvature in the model.

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

The predicted value for row nrow.

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

The derivative with respect to the j-th unknown parameter.

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

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

integer, intent(out) :: msg(nq,j)

The error checking results.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

The number of function evaluations.

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

A work array of (n, m, nq) elements.

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

A work array of (n, nq) elements.

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

A work array of (n, np, nq) elements.

public subroutine djckm(fcn, n, m, np, nq, beta, xplusd, ifixb, ifixx, ldifx, eta, tol, nrow, epsmac, j, lq, typj, h0, hc0, iswrtb, pv, d, diffj, msg1, msg, istop, nfev, wrk1, wrk2, wrk6, interval)

Check user supplied analytic derivatives against numerical derivatives.

Arguments

Type IntentOptional Attributes Name
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(inout) :: xplusd(n,m)

The values of x + delta.

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

The relative noise in the model.

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

The agreement tolerance.

integer, intent(in) :: nrow

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

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

The value of machine precision.

integer, intent(in) :: j

The index of the partial derivative being examined.

integer, intent(in) :: lq

The response currently being examined.

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

The typical size of the j-th unknown beta or delta.

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

The initial step size for the finite difference derivative.

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

The relative step size for central finite differences.

logical, intent(in) :: iswrtb

The variable designating whether the derivatives wrt beta (iswrtb = .true.) or delta (iswrtb = .false.) are being checked.

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

The predicted value for row nrow.

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

The derivative with respect to the j-th unknown parameter.

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

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

integer, intent(out) :: msg1

The first set of error checking results.

integer, intent(out) :: msg(nq,j)

The error checking results.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

The number of function evaluations.

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

A work array of (n, m, nq) elements.

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

A work array of (n, nq) elements.

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

A work array of (n, np, nq) elements.

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

Specifies which checks can be performed when checking derivatives based on the interval of the bound constraints.

public subroutine djckz(fcn, n, m, np, nq, beta, xplusd, ifixb, ifixx, ldifx, nrow, epsmac, j, lq, iswrtb, tol, d, fd, typj, pvpstp, stp0, pv, diffj, msg, istop, nfev, wrk1, wrk2, wrk6)

Recheck the derivatives in the case where the finite difference derivative disagrees with the analytic derivative and the analytic derivative is zero.

Arguments

Type IntentOptional Attributes Name
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(inout) :: xplusd(n,m)

The values of x + delta.

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

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

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

The value of machine precision.

integer, intent(in) :: j

The index of the partial derivative being examined.

integer, intent(in) :: lq

The response currently being examined.

logical, intent(in) :: iswrtb

The variable designating whether the derivatives wrt beta (iswrtb = .true.) or delta (iswrtb = .false.) are being checked.

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

The agreement tolerance.

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

The derivative with respect to the j-th unknown parameter.

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

The forward difference derivative wrt the j-th parameter.

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

The typical size of the j-th unknown beta or delta.

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

The predicted value for row nrow of the model using the current parameter estimates for all but the j-th parameter value, which is beta(j) + stp0.

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

The initial step size for the finite difference derivative.

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

The predicted value from the model for row nrow.

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

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

integer, intent(out) :: msg(nq,j)

The error checking results.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

The number of function evaluations.

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

A work array of (n, m, nq) elements.

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

A work array of (n, nq) elements.

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

A work array of (n, np, nq) elements.

public pure subroutine dodchk(n, m, np, nq, isodr, anajac, implct, beta, ifixb, ldx, ldifx, ldscld, ldstpd, ldwe, ld2we, ldwd, ld2wd, ldy, lwork, lwkmn, liwork, liwkmn, sclb, scld, stpb, stpd, info, lower, upper)

Check input parameters, indicating errors found using nonzero values of argument info.

Arguments

Type IntentOptional Attributes Name
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.

logical, intent(in) :: isodr

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

logical, intent(in) :: anajac

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

logical, intent(in) :: implct

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

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.

integer, intent(in) :: ldx

The leading dimension of array x.

integer, intent(in) :: ldifx

The leading dimension of array ifixx.

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

The leading dimension of array y.

integer, intent(in) :: lwork

The length of vector work.

integer, intent(in) :: lwkmn

The minimum acceptable length of array work.

integer, intent(in) :: liwork

The length of vector iwork.

integer, intent(in) :: liwkmn

The minimum acceptable length of array iwork.

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.

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

The step for the finite difference derivative wrt beta.

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

The step for the finite difference derivative wrt delta.

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 subroutine dodstp(n, m, np, nq, npp, f, fjacb, fjacd, wd, ldwd, ld2wd, ss, tt, ldtt, delta, alpha, epsfcn, isodr, tfjacb, omega, u, qraux, kpvt, s, t, phi, irank, rcond, forvcv, wrk1, wrk2, wrk3, wrk4, wrk5, wrk, lwrk, tempret, istopc)

Compute locally constrained steps s and t, and phi(alpha).

Arguments

Type IntentOptional Attributes Name
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.

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

The (weighted) estimated values of epsilon.

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) :: wd(ldwd,ld2wd,m)

The (squared) 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) :: ss(np)

The scaling values for the unfixed betas.

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

The estimated errors in the explanatory variables.

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

The Levenberg-Marquardt parameter.

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

The function's precision.

logical, intent(in) :: isodr

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

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

The array omega*fjacb.

real(kind=wp), intent(out) :: omega(nq,nq)

The array defined such that: omega*trans(omega) = inv(I + fjacd*inv(e)*trans(fjacd)) = (I - fjacd*inv(p)*trans(fjacd)) where e = d**2 + alpha*tt**2 and p = trans(fjacd)*fjacd + d**2 + alpha*tt**2.

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

The approximate null vector for tfjacb.

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

The array required to recover the orthogonal part of the Q-R decomposition.

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

The pivot vector.

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

The step for beta.

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

The step for delta.

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

The difference between the norm of the scaled step and the trust region diameter.

integer, intent(out) :: irank

The rank deficiency of the Jacobian wrt beta.

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

The approximate reciprocal condition number of tfjacb.

logical, intent(in) :: forvcv

The variable designating whether this subroutine was called to set up for the covariance matrix computations (forvcv = .true.) or not (forvcv = .false.).

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

A work array of (n, nq, m) elements.

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

A work array of (n, nq) elements.

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

A work array of (np) elements.

real(kind=wp), intent(out) :: wrk4(m,m)

A work array of (m, m) elements.

real(kind=wp), intent(out) :: wrk5(m)

A work array of (m) elements.

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

A work array of (lwrk) elements, equivalenced to wrk1 and wrk2.

integer, intent(in) :: lwrk

The length of vector wrk.

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

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

integer, intent(inout) :: istopc

The variable designating whether the computations were stopped due to a numerical error within subroutine dodstp.

public subroutine dodvcv(n, m, np, nq, npp, f, fjacb, fjacd, wd, ldwd, ld2wd, ssf, ss, tt, ldtt, delta, epsfcn, isodr, vcv, sd, wrk6, omega, u, qraux, jpvt, s, t, irank, rcond, rss, idf, rvar, ifixb, wrk1, wrk2, wrk3, wrk4, wrk5, wrk, lwrk, tempret, istopc)

Compute covariance matrix of estimated parameters.

Arguments

Type IntentOptional Attributes Name
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.

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

The (weighted) estimated values of epsilon.

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

The scaling values used for beta.

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

The scaling values for the unfixed betas.

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

The estimated errors in the explanatory variables.

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

The function's precision.

logical, intent(in) :: isodr

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

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

The covariance matrix of the estimated betas.

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

The standard deviations of the estimated betas.

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

A work array of (n*nq, np) elements.

real(kind=wp), intent(out) :: omega(nq,nq)

The array defined such that omega*trans(omega) = inv(I + fjacd*inv(e)*trans(fjacd)) = (I - fjacd*inv(p)*trans(fjacd)).

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

The approximate null vector for fjacb.

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

The array required to recover the orthogonal part of the Q-R decomposition.

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

The pivot vector.

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

The step for beta.

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

The step for delta.

integer, intent(out) :: irank

The rank deficiency of the Jacobian wrt beta.

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

The approximate reciprocal condition of fjacb.

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

The residual sum of squares.

integer, intent(out) :: 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(out) :: rvar

The residual variance.

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(out) :: wrk1(n,nq,m)

A work array of (n, nq, m) elements.

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

A work array of (n, nq) elements.

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

A work array of (np) elements.

real(kind=wp), intent(out) :: wrk4(m,m)

A work array of (m, m) elements.

real(kind=wp), intent(out) :: wrk5(m)

A work array of (m) elements.

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

A work array of (lwrk) elements, equivalenced to wrk1 and wrk2.

integer, intent(in) :: lwrk

The length of vector lwrk.

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

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

integer, intent(out) :: istopc

The variable designating whether the computations were stoped due to a numerical error within subroutine dodstp.

public pure subroutine dpack(n2, n1, v1, v2, ifix)

Select the unfixed elements of v2 and return them in v1.

Arguments

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

The number of items in v2.

integer, intent(out) :: n1

The number of items in v1.

real(kind=wp), intent(out) :: v1(n2)

The vector of the unfixed items from v2.

real(kind=wp), intent(in) :: v2(n2)

The vector of the fixed and unfixed items from which the unfixed elements are to be extracted.

integer, intent(in) :: ifix(n2)

The values designating whether the elements of v2 are fixed at their input values or not.

public subroutine dpvb(fcn, n, m, np, nq, beta, xplusd, ifixb, ifixx, ldifx, nrow, j, lq, stp, istop, nfev, pvb, wrk1, wrk2, wrk6)

Compute the nrow-th function value using beta(j) + stp.

Arguments

Type IntentOptional Attributes Name
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 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) :: xplusd(n,m)

The values of x + delta.

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

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

integer, intent(in) :: j

The index of the partial derivative being examined.

integer, intent(in) :: lq

The response currently being examined.

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

The step size for the finite difference derivative.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

The number of function evaluations.

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

The function value for the selected observation & response.

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

Work array.

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

Work array.

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

Work array.

public subroutine dpvd(fcn, n, m, np, nq, beta, xplusd, ifixb, ifixx, ldifx, nrow, j, lq, stp, istop, nfev, pvd, wrk1, wrk2, wrk6)

Compute nrow-th function value using x(nrow, j) + delta(nrow, j) + stp.

Arguments

Type IntentOptional Attributes Name
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 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(in) :: beta(np)

The function parameters.

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

The values of x + delta.

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

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

integer, intent(in) :: j

The index of the partial derivative being examined.

integer, intent(in) :: lq

The response currently being examined.

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

The step size for the finite difference derivative.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

The number of function evaluations.

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

The function value for the selected observation & response.

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

Work array.

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

Work array.

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

Work array.

public pure subroutine dscale(n, m, scl, ldscl, t, ldt, sclt, ldsclt)

Scale t by the inverse of scl, i.e., compute t/scl.

Arguments

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

The number of rows of data in t.

integer, intent(in) :: m

The number of columns of data in t.

real(kind=wp), intent(in) :: scl(ldscl,m)

The scale values.

integer, intent(in) :: ldscl

The leading dimension of array scl.

real(kind=wp), intent(in) :: t(ldt,m)

The array to be inversely scaled by scl.

integer, intent(in) :: ldt

The leading dimension of array t.

real(kind=wp), intent(out) :: sclt(ldsclt,m)

The inversely scaled matrix.

integer, intent(in) :: ldsclt

The leading dimension of array sclt.

public pure subroutine dsclb(np, beta, ssf)

Select scaling values for beta according to the algorithm given in the ODRPACK95 reference guide.

Arguments

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

The number of function parameters.

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

The function parameters.

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

The scaling values for beta.

public pure subroutine dscld(n, m, x, ldx, tt, ldtt)

Select scaling values for delta according to the algorithm given in the ODRPACK95 reference guide.

Arguments

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

The number of observations.

integer, intent(in) :: m

The number of columns of data in the independent variable.

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(out) :: tt(ldtt,m)

The scaling values for delta.

integer, intent(in) :: ldtt

The leading dimension of array tt.

public pure subroutine dsetn(n, m, x, ldx, nrow)

Select the row at which the derivative will be checked.

Arguments

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

The number of observations.

integer, intent(in) :: m

The number of columns of data in the independent variable.

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

The independent variable.

integer, intent(in) :: ldx

The leading dimension of array x.

integer, intent(inout) :: nrow

The selected row number of the independent variable.

public pure subroutine dsolve(n, t, ldt, b, job)

Solve systems of the form:

Read more…

Arguments

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

The number of rows and columns of data in array t.

real(kind=wp), intent(in) :: t(ldt,n)

The upper or lower tridiagonal system.

integer, intent(in) :: ldt

The leading dimension of array t.

real(kind=wp), intent(inout) :: b(n)

On input: the right hand side; On exit: the solution.

integer, intent(in) :: job

What kind of system is to be solved: 1: Solve t * x = b, where t is lower triangular, 2: Solve t * x = b, where t is upper triangular, 3: Solve trans(t) * x = b, where t is lower triangular, 4: Solve trans(t) * x = b, where t is upper triangular.

public pure subroutine dunpac(n2, v1, v2, ifix)

Copy the elements of v1 into the locations of v2 which are unfixed.

Arguments

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

The number of items in v2.

real(kind=wp), intent(in) :: v1(n2)

The vector of the unfixed items.

real(kind=wp), intent(out) :: v2(n2)

The vector of the fixed and unfixed items into which the elements of v1 are to be inserted.

integer, intent(in) :: ifix(n2)

The values designating whether the elements of v2 are fixed at their input values or not.

public pure subroutine dvevtr(m, nq, indx, v, ldv, ld2v, e, lde, ve, ldve, ld2ve, vev, ldvev, wrk5)

Compute v*e*trans(v) for the (indx)th m by nq array in v.

Arguments

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

The number of columns of data in the independent variable.

integer, intent(in) :: nq

The number of responses per observation.

integer, intent(in) :: indx

The row in v in which the m by nq array is stored.

real(kind=wp), intent(in) :: v(ldv,ld2v,nq)

An array of nq by m matrices.

integer, intent(in) :: ldv

The leading dimension of array v.

integer, intent(in) :: ld2v

The second dimension of array v.

real(kind=wp), intent(in) :: e(lde,m)

The m by m matrix of the factors, so ete = (d**2 + alpha*t**2).

integer, intent(in) :: lde

The leading dimension of array e.

real(kind=wp), intent(out) :: ve(ldve,ld2ve,m)

The nq by m array ve = v * inv(e).

integer, intent(in) :: ldve

The leading dimension of array ve.

integer, intent(in) :: ld2ve

The second dimension of array ve.

real(kind=wp), intent(out) :: vev(ldvev,nq)

The nq by nq array vev = v * inv(ete) * trans(v).

integer, intent(in) :: ldvev

The leading dimension of array vev.

real(kind=wp), intent(out) :: wrk5(m)

An m work vector.

public pure subroutine dwght(n, m, wt, ldwt, ld2wt, t, wtt)

Scale matrix t using wt, i.e., compute wtt = wt*t.

Arguments

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

The number of rows of data in t.

integer, intent(in) :: m

The number of columns of data in t.

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

The weights.

integer, intent(in) :: ldwt

The leading dimension of array wt.

integer, intent(in) :: ld2wt

The second dimension of array wt.

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

The array being scaled by wt.

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

The results of weighting array t by wt. Array wtt can be the same as t only if the arrays in wt are upper triangular with zeros below the diagonal.

public pure subroutine dwinf(n, m, np, nq, ldwe, ld2we, isodr, deltai, epsi, xplusi, fni, sdi, vcvi, rvari, wssi, wssdei, wssepi, rcondi, etai, olmavi, taui, alphai, actrsi, pnormi, rnorsi, prersi, partli, sstoli, taufci, epsmai, beta0i, betaci, betasi, betani, si, ssi, ssfi, qrauxi, ui, fsi, fjacbi, we1i, diffi, deltsi, deltni, ti, tti, omegai, fjacdi, wrk1i, wrk2i, wrk3i, wrk4i, wrk5i, wrk6i, wrk7i, loweri, upperi, lwkmn)

Get storage locations within REAL (wp) work space.

Arguments

Type IntentOptional Attributes Name
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) :: ldwe

The leading dimension of array we.

integer, intent(in) :: ld2we

The second dimension of array we.

logical, intent(in) :: isodr

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

integer, intent(out) :: deltai

The starting location in array work of array delta.

integer, intent(out) :: epsi

The starting location in array work of array eps.

integer, intent(out) :: xplusi

The starting location in array work of array xplusd.

integer, intent(out) :: fni

The starting location in array work of array fn.

integer, intent(out) :: sdi

The starting location in array work of array sd.

integer, intent(out) :: vcvi

The starting location in array work of array vcv.

integer, intent(out) :: rvari

The location in array work of variable rvar.

integer, intent(out) :: wssi

The location in array work of variable wss.

integer, intent(out) :: wssdei

The location in array work of variable wssdel.

integer, intent(out) :: wssepi

The location in array work of variable wsseps.

integer, intent(out) :: rcondi

The location in array work of variable rcondi.

integer, intent(out) :: etai

The location in array work of variable eta.

integer, intent(out) :: olmavi

The location in array work of variable olmavg.

integer, intent(out) :: taui

The location in array work of variable tau.

integer, intent(out) :: alphai

The location in array work of variable alpha.

integer, intent(out) :: actrsi

The location in array work of variable actrs.

integer, intent(out) :: pnormi

The location in array work of variable pnorm.

integer, intent(out) :: rnorsi

The location in array work of variable rnorms.

integer, intent(out) :: prersi

The location in array work of variable prers.

integer, intent(out) :: partli

The location in array work of variable partol.

integer, intent(out) :: sstoli

The location in array work of variable sstol.

integer, intent(out) :: taufci

The location in array work of variable taufac.

integer, intent(out) :: epsmai

The location in array work of variable epsmac.

integer, intent(out) :: beta0i

The starting location in array work of array beta0.

integer, intent(out) :: betaci

The starting location in array work of array betac.

integer, intent(out) :: betasi

The starting location in array work of array betas.

integer, intent(out) :: betani

The starting location in array work of array betan.

integer, intent(out) :: si

The starting location in array work of array s.

integer, intent(out) :: ssi

The starting location in array work of array ss.

integer, intent(out) :: ssfi

The starting location in array work of array ssf.

integer, intent(out) :: qrauxi

The starting location in array work of array qraux.

integer, intent(out) :: ui

The starting location in array work of array u.

integer, intent(out) :: fsi

The starting location in array work of array fs.

integer, intent(out) :: fjacbi

The starting location in array work of array fjacb.

integer, intent(out) :: we1i

The starting location in array work of array we1.

integer, intent(out) :: diffi

The starting location in array work of array diff.

integer, intent(out) :: deltsi

The starting location in array work of array deltas.

integer, intent(out) :: deltni

The starting location in array work of array deltan.

integer, intent(out) :: ti

The starting location in array work of array t.

integer, intent(out) :: tti

The starting location in array work of array tt.

integer, intent(out) :: omegai

The starting location in array work of array omega.

integer, intent(out) :: fjacdi

The starting location in array work of array fjacd.

integer, intent(out) :: wrk1i

The starting location in array work of array wrk1.

integer, intent(out) :: wrk2i

The starting location in array work of array wrk2.

integer, intent(out) :: wrk3i

The starting location in array work of array wrk3.

integer, intent(out) :: wrk4i

The starting location in array work of array wrk4.

integer, intent(out) :: wrk5i

The starting location in array work of array wrk5.

integer, intent(out) :: wrk6i

The starting location in array work of array wrk6.

integer, intent(out) :: wrk7i

The starting location in array work of array wrk7.

integer, intent(out) :: loweri

The starting location in array work of array lower.

integer, intent(out) :: upperi

The starting location in array work of array upper.

integer, intent(out) :: lwkmn

The minimum acceptable length of vector work.

public pure subroutine mbfb(np, beta, lower, upper, ssf, stpb, neta, eta, interval)

Ensure range of bounds is large enough for derivative checking. Move beta away from bounds so that derivatives can be calculated.

Arguments

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

The number of parameters np.

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

Function parameters.

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

!! Lower bound on beta.

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

Upper bound on beta.

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

The scale used for the betas.

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

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

integer, intent(in) :: neta

Number of good digits in the function results.

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

The relative noise in the function results.

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

Specifies which difference methods and step sizes are supported by the current interval upper-lower.