odr_long_c Subroutine

public subroutine odr_long_c(fcn, n, m, q, np, ldwe, ld2we, ldwd, ld2wd, ldifx, ldstpd, ldscld, lrwork, liwork, beta, y, x, we, wd, ifixb, ifixx, stpb, stpd, sclb, scld, delta, lower, upper, rwork, iwork, job, ndigit, taufac, sstol, partol, maxit, iprint, lunerr, lunrpt, info) bind(C)

"Long-call" wrapper for the odr routine including all optional arguments.

Arguments

Type IntentOptional Attributes Name
procedure(fcn_tc) :: fcn

User-supplied subroutine for evaluating the model.

integer(kind=c_int), intent(in) :: n

Number of observations.

integer(kind=c_int), intent(in) :: m

Number of columns of data in the independent variable.

integer(kind=c_int), intent(in) :: q

Number of responses per observation.

integer(kind=c_int), intent(in) :: np

Number of function parameters.

integer(kind=c_int), intent(in) :: ldwe

Leading dimension of array we, ldwe ∈ {1, n}.

integer(kind=c_int), intent(in) :: ld2we

Second dimension of array we, ld2we ∈ {1, q}.

integer(kind=c_int), intent(in) :: ldwd

Leading dimension of array wd, ldwd ∈ {1, n}.

integer(kind=c_int), intent(in) :: ld2wd

Second dimension of array wd, ld2wd ∈ {1, m}.

integer(kind=c_int), intent(in) :: ldifx

Leading dimension of array ifixx, ldifx ∈ {1, n}.

integer(kind=c_int), intent(in) :: ldstpd

Leading dimension of array stpd, ldstpd ∈ {1, n}.

integer(kind=c_int), intent(in) :: ldscld

Leading dimension of array scld, ldscld ∈ {1, n}.

integer(kind=c_int), intent(in) :: lrwork

Length of array rwork.

integer(kind=c_int), intent(in) :: liwork

Length of array iwork.

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

Function parameters.

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

Dependent variable. Unused when the model is implicit.

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

Explanatory variable.

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

epsilon weights.

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

delta weights.

integer(kind=c_int), intent(in), optional :: ifixb(np)

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

integer(kind=c_int), intent(in), optional :: ifixx(ldifx,m)

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

real(kind=c_double), intent(in), optional :: stpb(np)

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

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

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

real(kind=c_double), intent(in), optional :: sclb(np)

Scaling values for beta.

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

Scaling values for delta.

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

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

real(kind=c_double), intent(in), optional :: lower(np)

Lower bound on beta.

real(kind=c_double), intent(in), optional :: upper(np)

Upper bound on beta.

real(kind=c_double), intent(inout), optional :: rwork(lrwork)

Real work space.

integer(kind=c_int), intent(inout), optional :: iwork(liwork)

Integer work space.

integer(kind=c_int), intent(in), optional :: job

Variable controlling initialization and computational method.

integer(kind=c_int), intent(in), optional :: ndigit

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

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

Factor used to compute the initial trust region diameter.

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

Sum-of-squares convergence stopping tolerance.

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

Parameter convergence stopping tolerance.

integer(kind=c_int), intent(in), optional :: maxit

Maximum number of iterations allowed.

integer(kind=c_int), intent(in), optional :: iprint

Print control variable.

integer(kind=c_int), intent(in), optional :: lunerr

Logical unit number for error messages. 0: No output. 6: Output to standard output. k /= 0,6: Output to logical unit k.

integer(kind=c_int), intent(in), optional :: lunrpt

Logical unit number for computation reports. 0: No output. 6: Output to standard output. k /= 0,6: Output to logical unit k.

integer(kind=c_int), intent(out), optional :: info

Variable designating why the computations were stopped.


Calls

proc~~odr_long_c~~CallsGraph proc~odr_long_c odr_long_c proc~odr odr proc~odr_long_c->proc~odr interface~ddot ddot proc~odr->interface~ddot proc~access_workspace access_workspace proc~odr->proc~access_workspace proc~check_inputs check_inputs proc~odr->proc~check_inputs proc~check_jac check_jac proc~odr->proc~check_jac proc~derstep derstep proc~odr->proc~derstep proc~eta_fcn eta_fcn proc~odr->proc~eta_fcn proc~eval_jac eval_jac proc~odr->proc~eval_jac proc~fctrw fctrw proc~odr->proc~fctrw proc~init_work init_work proc~odr->proc~init_work proc~loc_iwork loc_iwork proc~odr->proc~loc_iwork proc~loc_rwork loc_rwork proc~odr->proc~loc_rwork proc~move_beta move_beta proc~odr->proc~move_beta proc~pack_vec pack_vec proc~odr->proc~pack_vec proc~print_error_inputs print_error_inputs proc~odr->proc~print_error_inputs proc~print_errors print_errors proc~odr->proc~print_errors proc~print_header print_header proc~odr->proc~print_header proc~print_reports print_reports proc~odr->proc~print_reports proc~scale_mat scale_mat proc~odr->proc~scale_mat proc~select_row select_row proc~odr->proc~select_row proc~set_flags set_flags proc~odr->proc~set_flags proc~trust_region trust_region proc~odr->proc~trust_region proc~unpack_vec unpack_vec proc~odr->proc~unpack_vec proc~vcv_beta vcv_beta proc~odr->proc~vcv_beta proc~workspace_dimensions workspace_dimensions proc~odr->proc~workspace_dimensions proc~access_workspace->proc~loc_iwork proc~access_workspace->proc~loc_rwork proc~check_jac_value check_jac_value proc~check_jac->proc~check_jac_value proc~hstep hstep proc~check_jac->proc~hstep proc~derstep->proc~hstep proc~eval_jac->proc~scale_mat proc~eval_jac->proc~unpack_vec proc~jac_cdiff jac_cdiff proc~eval_jac->proc~jac_cdiff proc~jac_fwdiff jac_fwdiff proc~eval_jac->proc~jac_fwdiff proc~set_ifix set_ifix proc~eval_jac->proc~set_ifix proc~fctr fctr proc~fctrw->proc~fctr proc~init_work->proc~set_flags proc~scale_beta scale_beta proc~init_work->proc~scale_beta proc~scale_delta scale_delta proc~init_work->proc~scale_delta proc~move_beta->proc~hstep proc~print_errors->proc~print_error_inputs proc~print_errors->proc~print_header proc~print_error_derivative print_error_derivative proc~print_errors->proc~print_error_derivative proc~print_error_fcn print_error_fcn proc~print_errors->proc~print_error_fcn proc~print_reports->proc~print_header proc~print_reports->proc~set_flags proc~print_report_final print_report_final proc~print_reports->proc~print_report_final proc~print_report_initial print_report_initial proc~print_reports->proc~print_report_initial proc~print_report_iter print_report_iter proc~print_reports->proc~print_report_iter proc~trust_region->proc~scale_mat proc~lcstep lcstep proc~trust_region->proc~lcstep proc~scale_vec scale_vec proc~trust_region->proc~scale_vec dpodi dpodi proc~vcv_beta->dpodi proc~vcv_beta->proc~lcstep proc~check_jac_curv check_jac_curv proc~check_jac_value->proc~check_jac_curv proc~check_jac_zero check_jac_zero proc~check_jac_value->proc~check_jac_zero proc~fpvb fpvb proc~check_jac_value->proc~fpvb proc~fpvd fpvd proc~check_jac_value->proc~fpvd proc~jac_cdiff->proc~derstep proc~jac_cdiff->proc~hstep proc~jac_fwdiff->proc~derstep proc~jac_fwdiff->proc~hstep proc~lcstep->proc~scale_mat proc~lcstep->proc~fctr dchex dchex proc~lcstep->dchex dqrdc dqrdc proc~lcstep->dqrdc dqrsl dqrsl proc~lcstep->dqrsl dtrco dtrco proc~lcstep->dtrco dtrsl dtrsl proc~lcstep->dtrsl interface~drot drot proc~lcstep->interface~drot interface~drotg drotg proc~lcstep->interface~drotg proc~esubi esubi proc~lcstep->proc~esubi proc~solve_trl solve_trl proc~lcstep->proc~solve_trl proc~vevtr vevtr proc~lcstep->proc~vevtr proc~ppf_tstudent ppf_tstudent proc~print_report_final->proc~ppf_tstudent proc~print_report_initial->proc~hstep proc~check_jac_curv->proc~fpvb proc~check_jac_curv->proc~fpvd proc~check_jac_fp check_jac_fp proc~check_jac_curv->proc~check_jac_fp proc~check_jac_zero->proc~fpvb proc~check_jac_zero->proc~fpvd proc~ppf_normal ppf_normal proc~ppf_tstudent->proc~ppf_normal proc~vevtr->proc~solve_trl proc~check_jac_fp->proc~fpvb proc~check_jac_fp->proc~fpvd

Subroutines

subroutine fcn_(beta, xplusd, ifixb, ifixx, ideval, f, fjacb, fjacd, istop)

Arguments

Type IntentOptional Attributes Name
real(kind=c_double), intent(in) :: beta(:)
real(kind=c_double), intent(in) :: xplusd(:,:)
integer(kind=c_int), intent(in) :: ifixb(:)
integer(kind=c_int), intent(in) :: ifixx(:,:)
integer(kind=c_int), intent(in) :: ideval
real(kind=c_double), intent(out) :: f(:,:)
real(kind=c_double), intent(out) :: fjacb(:,:,:)
real(kind=c_double), intent(out) :: fjacd(:,:,:)
integer(kind=c_int), intent(out) :: istop

Source Code

   subroutine odr_long_c( &
      fcn, &
      n, m, q, np, &
      ldwe, ld2we, &
      ldwd, ld2wd, &
      ldifx, &
      ldstpd, ldscld, &
      lrwork, liwork, &
      beta, y, x, &
      we, wd, &
      ifixb, ifixx, &
      stpb, stpd, &
      sclb, scld, &
      delta, &
      lower, upper, &
      rwork, iwork, &
      job, ndigit, taufac, &
      sstol, partol, maxit, &
      iprint, lunerr, lunrpt, &
      info) bind(C)
   !! "Long-call" wrapper for the `odr` routine including all optional arguments.

      procedure(fcn_tc) :: fcn
         !! User-supplied subroutine for evaluating the model.
      integer(c_int), intent(in) :: n
         !! Number of observations.
      integer(c_int), intent(in) :: m
         !! Number of columns of data in the independent variable.
      integer(c_int), intent(in) :: q
         !! Number of responses per observation.
      integer(c_int), intent(in) :: np
         !! Number of function parameters.
      integer(c_int), intent(in) :: ldwe
         !! Leading dimension of array `we`, `ldwe ∈ {1, n}`.
      integer(c_int), intent(in) :: ld2we
         !! Second dimension of array `we`, `ld2we ∈ {1, q}`.
      integer(c_int), intent(in) :: ldwd
         !! Leading dimension of array `wd`, `ldwd ∈ {1, n}`.
      integer(c_int), intent(in) :: ld2wd
         !! Second dimension of array `wd`, `ld2wd ∈ {1, m}`.
      integer(c_int), intent(in) :: ldifx
         !! Leading dimension of array `ifixx`, `ldifx ∈ {1, n}`.
      integer(c_int), intent(in) :: ldstpd
         !! Leading dimension of array `stpd`, `ldstpd ∈ {1, n}`.
      integer(c_int), intent(in) :: ldscld
         !! Leading dimension of array `scld`, `ldscld ∈ {1, n}`.
      integer(c_int), intent(in) :: lrwork
         !! Length of array `rwork`.
      integer(c_int), intent(in) :: liwork
         !! Length of array `iwork`.
      real(c_double), intent(inout) :: beta(np)
         !! Function parameters.
      real(c_double), intent(in) :: y(n, q)
         !! Dependent variable. Unused when the model is implicit.
      real(c_double), intent(in) :: x(n, m)
         !! Explanatory variable.
      real(c_double), intent(in), optional :: we(ldwe, ld2we, q)
         !! `epsilon` weights.
      real(c_double), intent(in), optional :: wd(ldwd, ld2wd, m)
         !! `delta` weights.
      integer(c_int), intent(in), optional :: ifixb(np)
         !! Values designating whether the elements of `beta` are fixed at their input values or not.
      integer(c_int), intent(in), optional :: ifixx(ldifx, m)
         !! Values designating whether the elements of `x` are fixed at their input values or not.
      real(c_double), intent(in), optional :: stpb(np)
         !! Relative step for computing finite difference derivatives with respect to `beta`.
      real(c_double), intent(in), optional :: stpd(ldstpd, m)
         !! Relative step for computing finite difference derivatives with respect to `delta`.
      real(c_double), intent(in), optional :: sclb(np)
         !! Scaling values for `beta`.
      real(c_double), intent(in), optional :: scld(ldscld, m)
         !! Scaling values for `delta`.
      real(c_double), intent(inout), optional :: delta(n, m)
         !! Error in the `x` data. `Shape: (n, m)`. Initial guess on input and estimated value
         !! on output.
      real(c_double), intent(in), optional :: lower(np)
         !! Lower bound on `beta`.
      real(c_double), intent(in), optional :: upper(np)
         !! Upper bound on `beta`.
      real(c_double), intent(inout), optional :: rwork(lrwork)
         !! Real work space.
      integer(c_int), intent(inout), optional :: iwork(liwork)
         !! Integer work space.
      integer(c_int), intent(in), optional :: job
         !! Variable controlling initialization and computational method.
      integer(c_int), intent(in), optional :: ndigit
         !! Number of accurate digits in the function results, as supplied by the user.
      real(c_double), intent(in), optional :: taufac
         !! Factor used to compute the initial trust region diameter.
      real(c_double), intent(in), optional :: sstol
         !! Sum-of-squares convergence stopping tolerance.
      real(c_double), intent(in), optional :: partol
         !! Parameter convergence stopping tolerance.
      integer(c_int), intent(in), optional :: maxit
         !! Maximum number of iterations allowed.
      integer(c_int), intent(in), optional :: iprint
         !! Print control variable.
      integer(c_int), intent(in), optional :: lunerr
         !! Logical unit number for error messages.
         !!  `0`: No output.
         !!  `6`: Output to standard output.
         !!  `k /= 0,6`: Output to logical unit `k`.
      integer(c_int), intent(in), optional :: lunrpt
         !! Logical unit number for computation reports.
         !!  `0`: No output.
         !!  `6`: Output to standard output.
         !!  `k /= 0,6`: Output to logical unit `k`.
      integer(c_int), intent(out), optional :: info
         !! Variable designating why the computations were stopped.

      call odr(fcn_, n, m, q, np, beta, y, x, &
               we=we, wd=wd, &
               ifixb=ifixb, ifixx=ifixx, &
               delta=delta, &
               lower=lower, upper=upper, &
               job=job, &
               ndigit=ndigit, taufac=taufac, sstol=sstol, partol=partol, maxit=maxit, &
               iprint=iprint, lunerr=lunerr, lunrpt=lunrpt, &
               stpb=stpb, stpd=stpd, &
               sclb=sclb, scld=scld, &
               info=info, &
               rwork=rwork, iwork=iwork)

   contains

      subroutine fcn_(beta, xplusd, ifixb, ifixx, ideval, f, fjacb, fjacd, istop)
         real(c_double), intent(in) :: beta(:)
         real(c_double), intent(in) :: xplusd(:, :)
         integer(c_int), intent(in) :: ifixb(:)
         integer(c_int), intent(in) :: ifixx(:, :)
         integer(c_int), intent(in) :: ideval
         real(c_double), intent(out) :: f(:, :)
         real(c_double), intent(out) :: fjacb(:, :, :)
         real(c_double), intent(out) :: fjacd(:, :, :)
         integer(c_int), intent(out) :: istop

         call fcn(n, m, q, np, ldifx, &
                  beta, xplusd, ifixb, ifixx, ideval, f, fjacb, fjacd, istop)

      end subroutine fcn_

   end subroutine odr_long_c