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~odr odr proc~odr->module~odrpack_core proc~print_report_final print_report_final proc~print_report_final->module~odrpack_core proc~print_report_initial print_report_initial proc~print_report_initial->module~odrpack_core proc~print_reports print_reports proc~print_reports->module~odrpack_core

Abstract Interfaces

abstract interface

  • public subroutine fcn_t(beta, xplusd, ifixb, ifixx, ideval, f, fjacb, fjacd, istop)

    User-supplied subroutine for evaluating the model.

    Arguments

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

    Current values of parameters.

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

    Current value of explanatory variable, i.e., x + delta. Shape: (n, m).

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

    Indicators for "fixing" parameters (beta).

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

    Indicators for "fixing" explanatory variable (x). Shape: (ldifx, m).

    integer, intent(in) :: ideval

    Indicator for selecting computation to be performed.

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

    Predicted function values. Shape: (n, q).

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

    Jacobian with respect to beta. Shape: (n, np, q).

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

    Jacobian with respect to errors delta. Shape: (n, m, q).

    integer, intent(out) :: istop

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


Functions

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

Compute step size for center and forward difference calculations.

Arguments

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

Finite difference method being used. 0: Forward finite differences. 1: Central finite differences.

integer, intent(in) :: k

Index into beta where betak resides.

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

k-th function parameter.

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

Scale used for the betas.

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

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 hstep(itype, neta, i, j, stp, ldstp) result(res)

Set relative step size for finite difference derivatives.

Arguments

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

Finite difference method being used. 0: Forward finite differences. 1: Central finite differences.

integer, intent(in) :: neta

Number of good digits in the function results.

integer, intent(in) :: i

Identifier for selecting user-supplied step sizes.

integer, intent(in) :: j

Identifier for selecting user-supplied step sizes.

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

Step size for the finite difference derivative.

integer, intent(in) :: ldstp

Leading dimension of array stp.

Return Value real(kind=wp)

public pure function ppf_normal(p) result(res)

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

Probability at which the percent point is to be evaluated; 0 < p < 1.

Return Value real(kind=wp)

public pure function ppf_tstudent(p, idf) result(res)

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

Probability at which the percent point is to be evaluated; 0 < p < 1.

integer, intent(in) :: idf

Degrees of freedom.

Return Value real(kind=wp)


Subroutines

public subroutine trust_region(n, m, np, q, 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, 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

Number of observations.

integer, intent(in) :: m

Number of columns of data in the explanatory variable.

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

integer, intent(in) :: npp

Number of function parameters being estimated.

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

Weighted estimated values of epsilon.

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

Jacobian with respect to beta.

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

Jacobian with respect to delta.

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

delta weights.

integer, intent(in) :: ldwd

Leading dimension of array wd.

integer, intent(in) :: ld2wd

Second dimension of array wd.

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

Scaling values used for the unfixed betas.

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

Scale used for the deltas.

integer, intent(in) :: ldtt

Leading dimension of array tt.

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

Estimated errors in the explanatory variables.

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

Current Levenberg-Marquardt parameter.

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

Trust region diameter.

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

Function's precision.

logical, intent(in) :: isodr

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

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

Array omega*fjacb.

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

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

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

Approximate null vector for tfjacb.

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

Array required to recover the orthogonal part of the QR decomposition.

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

Pivot vector.

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

Step for beta.

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

Step for delta.

integer, intent(out) :: nlms

Number of Levenberg-Marquardt steps taken.

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

Approximate reciprocal condition of tfjacb.

integer, intent(out) :: irank

Aank deficiency of the Jacobian wrt beta.

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

Work array of (n, q, m) elements.

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

Work array of (n, q) elements.

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

Work array of (np) elements.

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

Work array of (m, m) elements.

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

Work array of (m) elements.

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

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

integer, intent(in) :: lwrk

Length of vector wrk.

integer, intent(out) :: istopc

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

public pure subroutine access_workspace(n, m, np, q, ldwe, ld2we, rwork, lrwork, 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

Number of observations.

integer, intent(in) :: m

Number of columns of data in the explanatory variable.

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

integer, intent(in) :: ldwe

Leading dimension of array we.

integer, intent(in) :: ld2we

Second dimension of array we.

real(kind=wp), intent(inout) :: rwork(lrwork)

Real work space.

integer, intent(in) :: lrwork

Length of vector rwork.

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

Integer work space.

integer, intent(in) :: liwork

Length of vector iwork.

logical, intent(in) :: access

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

logical, intent(in) :: isodr

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

integer, intent(out) :: jpvt

Pivot vector.

integer, intent(out) :: omega

Starting location in array rwork of array omega(q**2).

integer, intent(out) :: u

Starting location in array rwork of array u(np).

integer, intent(out) :: qraux

Starting location in array rwork of array qraux(np).

integer, intent(out) :: sd

Starting location in array rwork of array sd(np).

integer, intent(out) :: vcv

Starting location in array rwork of array vcv(np**2).

integer, intent(out) :: wrk1

Starting location in array rwork of array wrk1(n, m, q).

integer, intent(out) :: wrk2

Starting location in array rwork of array wrk2(n, q).

integer, intent(out) :: wrk3

Starting location in array rwork of array wrk3(np).

integer, intent(out) :: wrk4

Starting location in array rwork of array wrk4(m, m).

integer, intent(out) :: wrk5

Starting location in array rwork of array wrk5(m).

integer, intent(out) :: wrk6

Starting location in array rwork of array wrk6(n, np, q).

integer, intent(out) :: nnzw

Number of nonzero weighted observations.

integer, intent(out) :: npp

Number of function parameters actually estimated.

integer, intent(out) :: job

Variable controlling problem initialization and computational method.

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

Parameter convergence stopping tolerance.

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

Sum-of-squares convergence stopping tolerance.

integer, intent(out) :: maxit

Maximum number of iterations allowed.

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

Factor used to compute the initial trust region diameter.

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

Relative noise in the function results.

integer, intent(out) :: neta

Number of accurate digits in the function results.

integer, intent(out) :: lunrpt

Logical unit number used for computation reports.

integer, intent(out) :: ipr1

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

integer, intent(out) :: ipr2

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

integer, intent(out) :: ipr2f

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

integer, intent(out) :: ipr3

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

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

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

Residual variance, i.e. the standard deviation squared.

integer, intent(inout) :: idf

Degrees of freedom of the fit, equal to the number of observations with nonzero weighted derivatives minus the number of parameters being estimated.

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

Trust region diameter.

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

Levenberg-Marquardt parameter.

integer, intent(inout) :: niter

Number of iterations taken.

integer, intent(inout) :: nfev

Number of function evaluations.

integer, intent(inout) :: njev

Number of Jacobian evaluations.

integer, intent(inout) :: int2

Number of internal doubling steps.

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

Average number of Levenberg-Marquardt steps per iteration.

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

Approximate reciprocal condition of fjacb.

integer, intent(inout) :: irank

Rank deficiency of the Jacobian wrt beta.

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

Saved actual relative reduction in the sum-of-squares.

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

Norm of the scaled estimated parameters.

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

Saved predicted relative reduction in the sum-of-squares.

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

Norm of the saved weighted epsilons and deltas.

integer, intent(inout) :: istop

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

public pure subroutine esubi(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

Number of observations.

integer, intent(in) :: m

Number of columns of data in the independent variable.

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

Squared delta weights.

integer, intent(in) :: ldwd

Leading dimension of array wd.

integer, intent(in) :: ld2wd

Second dimension of array wd.

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

Levenberg-Marquardt parameter.

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

Scaling values used for delta.

integer, intent(in) :: ldtt

Leading dimension of array tt.

integer, intent(in) :: i

Indexing variable.

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

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

public subroutine eta_fcn(fcn, n, m, np, q, xplusd, beta, epsmac, nrow, partmp, pv0, ifixb, ifixx, ldifx, istop, nfev, eta, neta, fjacd, f, fjacb, wrk7, info, lower, upper)

Compute noise and number of good digits in function results.

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

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

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

Values of x + delta.

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

Function parameters.

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

Value of machine precision.

integer, intent(in) :: nrow

Row number at which the derivative is to be checked.

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

Model parameters.

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

Original predicted values.

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

Number of function evaluations.

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

Noise in the model results.

integer, intent(out) :: neta

Number of accurate digits in the model results.

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

Jacobian wrt delta.

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

Function values.

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

Jacobian wrt beta.

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

Work array of (5, q) elements.

integer, intent(out) :: info

Variable indicating the status of the computation.

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

Lower bound of beta.

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

Upper bound of beta.

public subroutine eval_jac(fcn, anajac, cdjac, n, m, np, q, betac, beta, stpb, ifixb, ifixx, ldifx, x, 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

User-supplied subroutine for evaluating the model.

logical, intent(in) :: anajac

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

logical, intent(in) :: cdjac

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

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

Number of responses per observation.

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

Current estimated values of the unfixed betas.

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

Function parameters.

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

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

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

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

Independent variable.

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

Estimated values of delta.

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

Values of x + delta.

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

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

integer, intent(in) :: ldstpd

Leading dimension of array stpd.

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

Scale used for the betas.

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

Scaling values used for delta.

integer, intent(in) :: ldtt

Leading dimension of array tt.

integer, intent(in) :: neta

Number of accurate digits in the function results.

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

Predicted values of the function at the current point.

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

Step used for computing finite difference derivatives with respect to delta.

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

Work array of (n, m, q) elements.

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

Work array of (n, q) elements.

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

Work array of (np) elements.

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

Work array of (n, np, q) 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,q)

Jacobian with respect to beta.

logical, intent(in) :: isodr

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

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

Jacobian with respect to delta.

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

Square roots of the epsilon weights in array we.

integer, intent(in) :: ldwe

Leading dimension of arrays we and we1.

integer, intent(in) :: ld2we

Second dimension of arrays we and we1.

integer, intent(inout) :: njev

Number of Jacobian evaluations.

integer, intent(inout) :: nfev

Number of function evaluations.

integer, intent(out) :: istop

Variable designating that the user wishes the computations stopped.

integer, intent(out) :: info

Variable designating why the computations were stopped.

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

Lower bound of beta.

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

Upper bound of beta.

public pure subroutine fctr(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

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

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

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

Leading dimension of array a.

integer, intent(in) :: n

Number of rows and columns of data in array a.

integer, intent(out) :: info

Output flag. 0: Factorization was completed. k: Error condition. The leading minor of order k is not positive (semi)definite.

public pure subroutine fctrw(n, m, q, 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

Number of observations.

integer, intent(in) :: m

Number of columns of data in the explanatory variable.

integer, intent(in) :: q

Number of responses per observation.

integer, intent(in) :: npp

Number of function parameters being estimated.

logical, intent(in) :: isodr

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

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

Squared epsilon weights.

integer, intent(in) :: ldwe

Leading dimension of array we.

integer, intent(in) :: ld2we

Second dimension of array we.

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

Squared delta weights.

integer, intent(in) :: ldwd

Leading dimension of array wd.

integer, intent(in) :: ld2wd

Second dimension of array wd.

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

Work array of (q, q) elements.

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

Work array of (m, m) elements.

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

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

integer, intent(out) :: nnzw

Number of nonzero weighted observations.

integer, intent(out) :: info

Variable designating why the computations were stopped.

public pure subroutine set_flags(job, restart, initd, dovcv, redoj, anajac, cdjac, chkjac, isodr, implicit)

Set flags indicating conditions specified by job.

Arguments

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

Variable controlling problem initialization and computational method.

logical, intent(out) :: restart

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

logical, intent(out) :: initd

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

logical, intent(out) :: dovcv

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

logical, intent(out) :: redoj

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

logical, intent(out) :: anajac

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

logical, intent(out) :: cdjac

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

logical, intent(out) :: chkjac

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

logical, intent(out) :: isodr

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

logical, intent(out) :: implicit

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

public pure subroutine set_ifix(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

Number of rows of data in the array.

integer, intent(in) :: m

Number of columns of data in the array.

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

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

integer, intent(in) :: ldifix

Leading dimension of array ifix.

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

Array being set to zero according to the elements of ifix.

integer, intent(in) :: ldt

Leading dimension of array t.

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

Resulting array.

integer, intent(in) :: ldtfix

Leading dimension of array tfix.

public pure subroutine loc_iwork(m, q, np, msgbi, msgdi, ifix2i, istopi, nnzwi, nppi, idfi, jobi, iprinti, lunerri, lunrpti, nrowi, ntoli, netai, maxiti, niteri, nfevi, njevi, int2i, iranki, ldtti, boundi, liwkmin)

Get storage locations within integer work space.

Arguments

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

Number of columns of data in the independent variable.

integer, intent(in) :: q

Number of responses per observation.

integer, intent(in) :: np

Number of function parameters.

integer, intent(out) :: msgbi

Starting location in array iwork of array msgb.

integer, intent(out) :: msgdi

Starting location in array iwork of array msgd.

integer, intent(out) :: ifix2i

Starting location in array iwork of array ifix2.

integer, intent(out) :: istopi

Location in array iwork of variable istop.

integer, intent(out) :: nnzwi

Location in array iwork of variable nnzw.

integer, intent(out) :: nppi

Location in array iwork of variable npp.

integer, intent(out) :: idfi

Location in array iwork of variable idf.

integer, intent(out) :: jobi

Location in array iwork of variable job.

integer, intent(out) :: iprinti

Location in array iwork of variable iprint.

integer, intent(out) :: lunerri

Location in array iwork of variable lunerr.

integer, intent(out) :: lunrpti

Location in array iwork of variable lunrpt.

integer, intent(out) :: nrowi

Location in array iwork of variable nrow.

integer, intent(out) :: ntoli

Location in array iwork of variable ntol.

integer, intent(out) :: netai

Location in array iwork of variable neta.

integer, intent(out) :: maxiti

Location in array iwork of variable maxit.

integer, intent(out) :: niteri

Location in array iwork of variable niter.

integer, intent(out) :: nfevi

Location in array iwork of variable nfev.

integer, intent(out) :: njevi

Location in array iwork of variable njev.

integer, intent(out) :: int2i

Location in array iwork of variable int2.

integer, intent(out) :: iranki

Location in array iwork of variable irank.

integer, intent(out) :: ldtti

Location in array iwork of variable ldtt.

integer, intent(out) :: boundi

Location in array iwork of variable bound.

integer, intent(out) :: liwkmin

Minimum acceptable length of array iwork.

public pure subroutine loc_rwork(n, m, q, np, ldwe, ld2we, isodr, deltai, epsi, xplusdi, fni, sdi, vcvi, rvari, wssi, wssdeli, wssepsi, rcondi, etai, olmavgi, taui, alphai, actrsi, pnormi, rnormsi, prersi, partoli, sstoli, taufaci, epsmaci, beta0i, betaci, betasi, betani, si, ssi, ssfi, qrauxi, ui, fsi, fjacbi, we1i, diffi, deltasi, deltani, ti, tti, omegai, fjacdi, wrk1i, wrk2i, wrk3i, wrk4i, wrk5i, wrk6i, wrk7i, loweri, upperi, lrwkmin)

Get storage locations within real work space.

Arguments

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

Number of observations.

integer, intent(in) :: m

Number of columns of data in the explanatory variable.

integer, intent(in) :: q

Number of responses per observation.

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: ldwe

Leading dimension of array we.

integer, intent(in) :: ld2we

Second dimension of array we.

logical, intent(in) :: isodr

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

integer, intent(out) :: deltai

Starting location in array rwork of array delta(n, m).

integer, intent(out) :: epsi

Starting location in array rwork of array eps(n, q).

integer, intent(out) :: xplusdi

Starting location in array rwork of array xplusd(n, m).

integer, intent(out) :: fni

Starting location in array rwork of array fn(n, q).

integer, intent(out) :: sdi

Starting location in array rwork of array sd(np).

integer, intent(out) :: vcvi

Starting location in array rwork of array vcv(np**2).

integer, intent(out) :: rvari

Location in array rwork of variable rvar.

integer, intent(out) :: wssi

Location in array rwork of variable wss.

integer, intent(out) :: wssdeli

Location in array rwork of variable wssdel.

integer, intent(out) :: wssepsi

Location in array rwork of variable wsseps.

integer, intent(out) :: rcondi

Location in array rwork of variable rcondi.

integer, intent(out) :: etai

Location in array rwork of variable eta.

integer, intent(out) :: olmavgi

Location in array rwork of variable olmavg.

integer, intent(out) :: taui

Location in array rwork of variable tau.

integer, intent(out) :: alphai

Location in array rwork of variable alpha.

integer, intent(out) :: actrsi

Location in array rwork of variable actrs.

integer, intent(out) :: pnormi

Location in array rwork of variable pnorm.

integer, intent(out) :: rnormsi

Location in array rwork of variable rnorms.

integer, intent(out) :: prersi

Location in array rwork of variable prers.

integer, intent(out) :: partoli

Location in array rwork of variable partol.

integer, intent(out) :: sstoli

Location in array rwork of variable sstol.

integer, intent(out) :: taufaci

Location in array rwork of variable taufac.

integer, intent(out) :: epsmaci

Location in array rwork of variable epsmac.

integer, intent(out) :: beta0i

Starting location in array rwork of array beta0(np).

integer, intent(out) :: betaci

Starting location in array rwork of array betac(np).

integer, intent(out) :: betasi

Starting location in array rwork of array betas(np).

integer, intent(out) :: betani

Starting location in array rwork of array betan(np).

integer, intent(out) :: si

Starting location in array rwork of array s(np).

integer, intent(out) :: ssi

Starting location in array rwork of array ss(np).

integer, intent(out) :: ssfi

Starting location in array rwork of array ssf(np).

integer, intent(out) :: qrauxi

Starting location in array rwork of array qraux(np).

integer, intent(out) :: ui

Starting location in array rwork of array u(np).

integer, intent(out) :: fsi

Starting location in array rwork of array fs(n, q).

integer, intent(out) :: fjacbi

Starting location in array rwork of array fjacb(n, np, q).

integer, intent(out) :: we1i

Starting location in array rwork of array we1(ldwe, ldwe2, q).

integer, intent(out) :: diffi

Starting location in array rwork of array diff(q*(np + m)).

integer, intent(out) :: deltasi

Starting location in array rwork of array deltas(n, m).

integer, intent(out) :: deltani

Starting location in array rwork of array deltan(n, m).

integer, intent(out) :: ti

Starting location in array rwork of array t(n, m).

integer, intent(out) :: tti

Starting location in array rwork of array tt(n, m).

integer, intent(out) :: omegai

Starting location in array rwork of array omega(q**2).

integer, intent(out) :: fjacdi

Starting location in array rwork of array fjacd(n, m, q).

integer, intent(out) :: wrk1i

Starting location in array rwork of array wrk1(n, m, q).

integer, intent(out) :: wrk2i

Starting location in array rwork of array wrk2(n, q).

integer, intent(out) :: wrk3i

Starting location in array rwork of array wrk3(np).

integer, intent(out) :: wrk4i

Starting location in array rwork of array wrk4(m, m).

integer, intent(out) :: wrk5i

Starting location in array rwork of array wrk5(m).

integer, intent(out) :: wrk6i

Starting location in array rwork of array wrk6(n, np, q).

integer, intent(out) :: wrk7i

Starting location in array rwork of array wrk7(5*q).

integer, intent(out) :: loweri

Starting location in array rwork of array lower(np).

integer, intent(out) :: upperi

Starting location in array rwork of array upper(np).

integer, intent(out) :: lrwkmin

Minimum acceptable length of vector rwork.

public pure subroutine init_work(n, m, np, rwork, lrwork, iwork, liwork, x, 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

Number of observations.

integer, intent(in) :: m

Number of columns of data in the independent variable.

integer, intent(in) :: np

Number of function parameters.

real(kind=wp), intent(out) :: rwork(lrwork)

Real work space.

integer, intent(in) :: lrwork

Length of vector rwork.

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

Integer work space.

integer, intent(in) :: liwork

Length of vector iwork.

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

Independent variable.

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

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

Scaling values for delta.

integer, intent(in) :: ldscld

Leading dimension of array scld.

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

Function parameters.

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

Scaling values for beta.

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

Sum-of-squares convergence stopping criteria.

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

Parameter convergence stopping criteria.

integer, intent(in) :: maxit

Maximum number of iterations allowed.

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

Factor used to compute the initial trust region diameter.

integer, intent(in) :: job

Variable controlling problem initialization and computational method.

integer, intent(in) :: iprint

Print control variable.

integer, intent(in) :: lunerr

Logical unit number used for error messages.

integer, intent(in) :: lunrpt

Logical unit number used for computation reports.

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

Lower bounds for the function parameters.

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

Upper bounds for the function parameters.

integer, intent(in) :: epsmai

Location in array rwork of variable epsmac.

integer, intent(in) :: sstoli

Location in array rwork of variable sstol.

integer, intent(in) :: partli

Location in array rwork of variable partol.

integer, intent(in) :: maxiti

Location in array iwork of variable maxit.

integer, intent(in) :: taufci

Location in array rwork of variable taufac.

integer, intent(in) :: jobi

Location in array iwork of variable job.

integer, intent(in) :: iprini

Location in array iwork of variable iprint.

integer, intent(in) :: luneri

Location in array iwork of variable lunerr.

integer, intent(in) :: lunrpi

Location in array iwork of variable lunrpt.

integer, intent(in) :: ssfi

Starting location in array rwork of array ssf.

integer, intent(in) :: tti

Starting location in array rwork of the array tt.

integer, intent(in) :: ldtti

Leading dimension of array tt.

integer, intent(in) :: deltai

Starting location in array rwork of array delta.

integer, intent(in) :: loweri

Starting location in array iwork of array lower.

integer, intent(in) :: upperi

Starting location in array iwork of array upper.

integer, intent(in) :: boundi

Location in array iwork of variable bound.

public subroutine jac_cdiff(fcn, n, m, np, q, beta, x, 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

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

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

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

Function parameters.

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

Explanatory variable.

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

Estimated errors in the explanatory variables.

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

Values of x + delta.

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

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

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

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

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

integer, intent(in) :: ldstpd

Leading dimension of array stpd.

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

Scaling values used for beta.

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

Scaling values used for delta.

integer, intent(in) :: ldtt

Leading dimension of array tt.

integer, intent(in) :: neta

Number of good digits in the function results.

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

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

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

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

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

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

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

A work array of (n, q) elements.

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

A work array of (np) elements.

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

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

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

Jacobian with respect to beta.

logical, intent(in) :: isodr

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

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

Jacobian with respect to delta.

integer, intent(inout) :: nfev

Number of function evaluations.

integer, intent(out) :: istop

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

integer, intent(out) :: info

Variable designating why the computations were stopped.

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

Lower bound on beta.

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

Upper bound on beta.

public subroutine jac_fwdiff(fcn, n, m, np, q, beta, x, 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

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

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

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

Function parameters.

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

Explanatory variable.

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

Estimated errors in the explanatory variables.

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

Values of x + delta.

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

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

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

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

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

integer, intent(in) :: ldstpd

Leading dimension of array stpd.

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

Scaling values used for beta.

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

Scaling values used for delta.

integer, intent(in) :: ldtt

Leading dimension of array tt.

integer, intent(in) :: neta

Number of good digits in the function results.

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

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

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

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

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

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

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

A work array of (n, q) elements.

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

A work array of (np) elements.

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

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

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

Jacobian with respect to beta.

logical, intent(in) :: isodr

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

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

Jacobian with respect to delta.

integer, intent(inout) :: nfev

Number of function evaluations.

integer, intent(out) :: istop

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

integer, intent(out) :: info

Variable designating why the computations were stopped.

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

Lower bound on beta.

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

Upper bound on beta.

public subroutine check_jac(fcn, n, m, np, q, 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

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

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

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

Function parameters.

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

Function parameters offset such that steps don't cross bounds.

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

Values of x + delta.

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

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

Step size for finite difference derivatives wrt beta.

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

Step size for finite difference derivatives wrt delta.

integer, intent(in) :: ldstpd

Leading dimension of array stpd.

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

Scaling values used for beta.

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

Scaling values used for delta.

integer, intent(in) :: ldtt

Leading dimension of array tt.

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

Relative noise in the function results.

integer, intent(in) :: neta

Number of reliable digits in the model results.

integer, intent(out) :: ntol

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

integer, intent(in) :: nrow

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

logical, intent(in) :: isodr

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

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

Value of machine precision.

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

Predicted values using the user supplied parameter estimates.

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

Jacobian with respect to beta.

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

Jacobian with respect to delta.

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

Error checking results for the Jacobian wrt beta.

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

Error checking results for the Jacobian wrt delta.

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

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

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

Number of function evaluations.

integer, intent(inout) :: njev

Number of Jacobian evaluations.

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

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

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

A work array of (n, q) elements.

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

A work array of (n, np, q) 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 check_jac_curv(fcn, n, m, np, q, 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

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

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

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

Function parameters.

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

Values of x + delta.

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

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

Relative noise in the model.

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

Agreement tolerance.

integer, intent(in) :: nrow

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

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

Value of machine precision.

integer, intent(in) :: j

Index of the partial derivative being examined.

integer, intent(in) :: lq

Response currently being examined.

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

Relative step size for central finite differences.

logical, intent(in) :: iswrtb

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

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

Forward difference derivative wrt the j-th parameter.

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

Typical size of the j-th unknown beta or delta.

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

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

Initial step size for the finite difference derivative.

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

Predicted value of the model for row nrow.

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

Derivative with respect to the j-th unknown parameter.

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

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

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

Error checking results.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

Number of function evaluations.

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

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

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

A work array of (n, q) elements.

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

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

public subroutine check_jac_fp(fcn, n, m, np, q, 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

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

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

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

Function parameters.

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

Values of x + delta.

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

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

Relative noise in the model.

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

Agreement tolerance.

integer, intent(in) :: nrow

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

integer, intent(in) :: j

Index of the partial derivative being examined.

integer, intent(in) :: lq

Response currently being examined.

logical, intent(in) :: iswrtb

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

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

Forward difference derivative wrt the j-th parameter.

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

Typical size of the j-th unknown beta or delta.

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

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

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

Predicted value for row nrow.

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

Derivative with respect to the j-th unknown parameter.

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

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

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

Error checking results.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

Number of function evaluations.

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

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

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

A work array of (n, q) elements.

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

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

public subroutine check_jac_value(fcn, n, m, np, q, 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

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

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

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

Function parameters.

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

Values of x + delta.

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

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

Relative noise in the model.

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

Agreement tolerance.

integer, intent(in) :: nrow

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

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

Value of machine precision.

integer, intent(in) :: j

Index of the partial derivative being examined.

integer, intent(in) :: lq

Response currently being examined.

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

Typical size of the j-th unknown beta or delta.

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

Initial step size for the finite difference derivative.

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

Relative step size for central finite differences.

logical, intent(in) :: iswrtb

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

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

Predicted value for row nrow.

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

Derivative with respect to the j-th unknown parameter.

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

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

integer, intent(out) :: msg1

First set of error checking results.

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

Error checking results.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

Number of function evaluations.

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

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

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

A work array of (n, q) elements.

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

A work array of (n, np, q) 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 check_jac_zero(fcn, n, m, np, q, 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

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

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

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

Function parameters.

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

Values of x + delta.

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

integer, intent(in) :: nrow

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

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

Value of machine precision.

integer, intent(in) :: j

Index of the partial derivative being examined.

integer, intent(in) :: lq

Response currently being examined.

logical, intent(in) :: iswrtb

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

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

Agreement tolerance.

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

Derivative with respect to the j-th unknown parameter.

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

Forward difference derivative wrt the j-th parameter.

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

Typical size of the j-th unknown beta or delta.

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

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

Initial step size for the finite difference derivative.

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

Predicted value from the model for row nrow.

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

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

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

Error checking results.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

Number of function evaluations.

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

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

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

A work array of (n, q) elements.

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

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

public pure subroutine check_inputs(n, m, np, q, isodr, anajac, beta, ifixb, ldifx, ldscld, ldstpd, ldwe, ld2we, ldwd, ld2wd, lrwork, lrwkmin, liwork, liwkmin, 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

Number of observations.

integer, intent(in) :: m

Number of columns of data in the explanatory variable.

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

logical, intent(in) :: isodr

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

logical, intent(in) :: anajac

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

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

Function parameters.

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

integer, intent(in) :: ldscld

Leading dimension of array scld.

integer, intent(in) :: ldstpd

Leading dimension of array stpd.

integer, intent(in) :: ldwe

Leading dimension of array we.

integer, intent(in) :: ld2we

Second dimension of array we.

integer, intent(in) :: ldwd

Leading dimension of array wd.

integer, intent(in) :: ld2wd

Second dimension of array wd.

integer, intent(in) :: lrwork

Length of vector rwork.

integer, intent(in) :: lrwkmin

Minimum acceptable length of array rwork.

integer, intent(in) :: liwork

Length of vector iwork.

integer, intent(in) :: liwkmin

Minimum acceptable length of array iwork.

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

Scaling values for beta.

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

Scaling value for delta.

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

Step for the finite difference derivative wrt beta.

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

Step for the finite difference derivative wrt delta.

integer, intent(out) :: info

Variable designating why the computations were stopped.

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

Lower bound on beta.

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

Upper bound on beta.

public subroutine lcstep(n, m, np, q, 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, istopc)

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

Arguments

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

Number of observations.

integer, intent(in) :: m

Number of columns of data in the explanatory variable.

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

integer, intent(in) :: npp

Number of function parameters being estimated.

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

Weighted estimated values of epsilon.

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

Jacobian with respect to beta.

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

Jacobian with respect to delta.

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

Squared delta weights.

integer, intent(in) :: ldwd

Leading dimension of array wd.

integer, intent(in) :: ld2wd

Second dimension of array wd.

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

Scaling values for the unfixed betas.

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

Scaling values for delta.

integer, intent(in) :: ldtt

Leading dimension of array tt.

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

Estimated errors in the explanatory variables.

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

Levenberg-Marquardt parameter.

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

Function's precision.

logical, intent(in) :: isodr

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

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

Array omega*fjacb.

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

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)

Approximate null vector for tfjacb.

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

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

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

Pivot vector.

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

Step for beta.

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

Step for delta.

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

Difference between the norm of the scaled step and the trust region diameter.

integer, intent(out) :: irank

Rank deficiency of the Jacobian wrt beta.

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

Approximate reciprocal condition number of tfjacb.

logical, intent(in) :: forvcv

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,q,m)

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

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

A work array of (n, q) 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

Length of vector wrk.

integer, intent(inout) :: istopc

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

public subroutine vcv_beta(n, m, np, q, 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, istopc)

Compute covariance matrix of estimated parameters.

Arguments

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

Number of observations.

integer, intent(in) :: m

Number of columns of data in the explanatory variable.

integer, intent(in) :: np

Number of function parameters.

integer, intent(in) :: q

Number of responses per observation.

integer, intent(in) :: npp

Number of function parameters being estimated.

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

Weighted estimated values of epsilon.

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

Jacobian with respect to beta.

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

Jacobian with respect to delta.

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

delta weights.

integer, intent(in) :: ldwd

Leading dimension of array wd.

integer, intent(in) :: ld2wd

Second dimension of array wd.

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

Scaling values used for beta.

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

Scaling values for the unfixed betas.

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

Scaling values for delta.

integer, intent(in) :: ldtt

Leading dimension of array tt.

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

Estimated errors in the explanatory variables.

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

Function's precision.

logical, intent(in) :: isodr

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

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

Covariance matrix of the estimated betas.

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

Standard deviations of the estimated betas.

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

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

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

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)

Approximate null vector for fjacb.

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

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

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

Pivot vector.

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

Step for beta.

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

Step for delta.

integer, intent(out) :: irank

Rank deficiency of the Jacobian wrt beta.

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

Approximate reciprocal condition of fjacb.

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

Residual sum of squares.

integer, intent(out) :: idf

Degrees of freedom of the fit, equal to the number of observations with nonzero weighted derivatives minus the number of parameters being estimated.

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

Residual variance.

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

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

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

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

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

A work array of (n, q) 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

Length of vector lwrk.

integer, intent(out) :: istopc

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

public pure subroutine pack_vec(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

Number of items in v2.

integer, intent(out) :: n1

Number of items in v1.

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

Vector of the unfixed items from v2.

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

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

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

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

public pure subroutine unpack_vec(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

Number of items in v2.

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

Vector of the unfixed items.

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

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

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

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

public subroutine fpvb(fcn, n, m, np, q, beta, xplusd, ifixb, ifixx, ldifx, nrow, j, lq, stp, istop, nfev, pvb, fjacd, f, fjacb)

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

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

Number of responses per observation.

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

Function parameters.

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

Values of x + delta.

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

integer, intent(in) :: nrow

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

integer, intent(in) :: j

Index of the partial derivative being examined.

integer, intent(in) :: lq

Response currently being examined.

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

Step size for the finite difference derivative.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

Number of function evaluations.

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

Function value for the selected observation & response.

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

Jacobian wrt delta.

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

Predicted function values.

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

Jocobian wrt beta.

public subroutine fpvd(fcn, n, m, np, q, beta, xplusd, ifixb, ifixx, ldifx, nrow, j, lq, stp, istop, nfev, pvd, fjacd, f, fjacb)

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

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

Number of responses per observation.

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

Function parameters.

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

Values of x + delta.

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

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

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

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

integer, intent(in) :: ldifx

Leading dimension of array ifixx.

integer, intent(in) :: nrow

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

integer, intent(in) :: j

Index of the partial derivative being examined.

integer, intent(in) :: lq

Response currently being examined.

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

Step size for the finite difference derivative.

integer, intent(out) :: istop

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

integer, intent(inout) :: nfev

Number of function evaluations.

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

Function value for the selected observation & response.

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

Jacobian wrt delta.

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

Predicted function values.

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

Jocobian wrt beta.

public pure subroutine scale_vec(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

Number of rows of data in t.

integer, intent(in) :: m

Number of columns of data in t.

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

Scale values.

integer, intent(in) :: ldscl

Leading dimension of array scl.

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

Array to be inversely scaled by scl.

integer, intent(in) :: ldt

Leading dimension of array t.

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

Inversely scaled matrix.

integer, intent(in) :: ldsclt

Leading dimension of array sclt.

public pure subroutine scale_beta(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

Number of function parameters.

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

Function parameters.

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

Scaling values for beta.

public pure subroutine scale_delta(n, m, x, 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

Number of observations.

integer, intent(in) :: m

Number of columns of data in the independent variable.

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

Independent variable.

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

Scaling values for delta.

integer, intent(in) :: ldtt

Leading dimension of array tt.

public pure subroutine select_row(n, m, x, nrow)

Select the row at which the derivative will be checked.

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.

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

Independent variable.

integer, intent(inout) :: nrow

Selected row number of the independent variable.

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

Solve systems of the form:

Read more…

Arguments

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

Number of rows and columns of data in array t.

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

Upper or lower tridiagonal system.

integer, intent(in) :: ldt

Leading dimension of array t.

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

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 vevtr(m, q, indx, v, ldv, ld2v, e, lde, ve, ldve, ld2ve, vev, ldvev, wrk5)

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

Arguments

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

Number of columns of data in the independent variable.

integer, intent(in) :: q

Number of responses per observation.

integer, intent(in) :: indx

Row in v in which the m by q array is stored.

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

An array of q by m matrices.

integer, intent(in) :: ldv

Leading dimension of array v.

integer, intent(in) :: ld2v

Second dimension of array v.

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

Matrix of the factors, so ete = (d**2 + alpha*t**2).

integer, intent(in) :: lde

Leading dimension of array e.

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

Array ve = v * inv(e).

integer, intent(in) :: ldve

Leading dimension of array ve.

integer, intent(in) :: ld2ve

Second dimension of array ve.

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

Array vev = v * inv(ete) * trans(v).

integer, intent(in) :: ldvev

Leading dimension of array vev.

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

Work vector.

public pure subroutine scale_mat(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

Number of rows of data in t.

integer, intent(in) :: m

Number of columns of data in t.

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

Array of shape conformable to (ldwt,ld2wt,m) holding the weights.

integer, intent(in) :: ldwt

Leading dimension of array wt.

integer, intent(in) :: ld2wt

Second dimension of array wt.

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

Array of shape conformable to (n,m) being scaled by wt.

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

Array of shape conformable to (n,m) holding the result of weighting array t by array 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 move_beta(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

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)

Scale used for the betas.

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

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

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.