"Short-call" wrapper for the odr
routine including very few optional arguments.
Type | Intent | Optional | 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. |
||
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(inout), | optional | :: | delta(n,m) |
Error in the |
|
real(kind=c_double), | intent(in), | optional | :: | lower(np) |
Lower bound on |
|
real(kind=c_double), | intent(in), | optional | :: | upper(np) |
Upper bound on |
|
integer(kind=c_int), | intent(in), | optional | :: | job |
Variable controlling initialization and computational method. |
Type | Intent | Optional | 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 |
subroutine odr_short_c( & fcn, & n, m, q, np, & beta, y, x, & delta, & lower, upper, & job) bind(C) !! "Short-call" wrapper for the `odr` routine including very few 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. 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(inout), optional :: delta(n, m) !! Error in the `x` data. 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`. integer(c_int), intent(in), optional :: job !! Variable controlling initialization and computational method. call odr(fcn_, n, m, q, np, beta, y, x, & delta=delta, & lower=lower, upper=upper, & job=job) 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, 1, & beta, xplusd, ifixb, ifixx, ideval, f, fjacb, fjacd, istop) end subroutine fcn_ end subroutine odr_short_c