Check input parameters, indicating errors found using nonzero values of argument info
as described in the ODRPACK95 reference guide.
Type | Intent | Optional | 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 ( |
||
real(kind=wp), | intent(in) | :: | we(ldwe,ld2we,q) |
Squared |
||
integer, | intent(in) | :: | ldwe |
Leading dimension of array |
||
integer, | intent(in) | :: | ld2we |
Second dimension of array |
||
real(kind=wp), | intent(in) | :: | wd(ldwd,ld2wd,m) |
Squared |
||
integer, | intent(in) | :: | ldwd |
Leading dimension of array |
||
integer, | intent(in) | :: | ld2wd |
Second dimension of array |
||
real(kind=wp), | intent(out) | :: | wrk0(q,q) |
Work array of |
||
real(kind=wp), | intent(out) | :: | wrk4(m,m) |
Work array of |
||
real(kind=wp), | intent(out) | :: | we1(ldwe,ld2we,q) |
Factored |
||
integer, | intent(out) | :: | nnzw |
Number of nonzero weighted observations. |
||
integer, | intent(out) | :: | info |
Variable designating why the computations were stopped. |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public | :: | i | ||||
integer, | public | :: | finfo | ||||
integer, | public | :: | j | ||||
logical, | public | :: | notzero | ||||
logical, | public | :: | exited |
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. use odrpack_kinds, only: zero 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(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(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(wp), intent(out) :: wrk0(q, q) !! Work array of `(q, q)` elements. real(wp), intent(out) :: wrk4(m, m) !! Work array of `(m, m)` elements. real(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. ! Local scalars integer :: i, finfo, j logical :: notzero, exited ! Variable Definitions (alphabetically) ! I: An indexing variable. ! J: An indexing variable. ! L: An indexing variable. ! NOTZERO: The variable designating whether a given component of the weight array WE ! contains a nonzero element (FALSE) or not (TRUE). ! Check EPSILON weights, and store factorization in WE1 exited = .false. if (we(1, 1, 1) < zero) then ! WE contains a scalar we1(1, 1, 1) = -sqrt(abs(we(1, 1, 1))) nnzw = n else nnzw = 0 if (ldwe == 1) then if (ld2we == 1) then ! WE contains a diagonal matrix do j = 1, q if (we(1, 1, j) > zero) then nnzw = n we1(1, 1, j) = sqrt(we(1, 1, j)) elseif (we(1, 1, j) < zero) then info = 30010 exited = .true. exit end if end do else ! WE contains a full Q by Q semidefinite matrix do j = 1, q wrk0(1:j, j) = we(1, 1:j, j) end do call fctr(.true., wrk0, q, q, finfo) if (finfo /= 0) then info = 30010 exited = .true. else do j = 1, q we1(1, :, j) = wrk0(:, j) if (we1(1, j, j) /= zero) then nnzw = n end if end do end if end if else if (ld2we == 1) then ! WE contains an array of diagonal matrix do i = 1, n notzero = .false. do j = 1, q if (we(i, 1, j) > zero) then notzero = .true. we1(i, 1, j) = sqrt(we(i, 1, j)) elseif (we(i, 1, j) < zero) then info = 30010 exited = .true. exit end if end do if (exited) exit if (notzero) then nnzw = nnzw + 1 end if end do else ! WE contains an array of full Q by Q semidefinite matrices do i = 1, n do j = 1, q wrk0(1:j, j) = we(i, 1:j, j) end do call fctr(.true., wrk0, q, q, finfo) if (finfo /= 0) then info = 30010 exited = .true. exit else notzero = .false. do j = 1, q we1(i, :, j) = wrk0(:, j) if (we1(i, j, j) /= zero) then notzero = .true. end if end do end if if (notzero) then nnzw = nnzw + 1 end if end do end if end if end if ! Check for a sufficient number of nonzero EPSILON weights if (.not. exited) then if (nnzw < npp) info = 30020 end if ! Check DELTA weights if (.not. isodr .or. wd(1, 1, 1) < zero) then ! Problem is not ODR, or WD contains a scalar return else if (ldwd == 1) then if (ld2wd == 1) then ! WD contains a diagonal matrix if (any(wd(1, 1, :) <= zero)) then info = max(30001, info + 1) return end if else ! WD contains a full M by M positive definite matrix do j = 1, m wrk4(1:j, j) = wd(1, 1:j, j) end do call fctr(.false., wrk4, m, m, finfo) if (finfo /= 0) then info = max(30001, info + 1) return end if end if else if (ld2wd == 1) then ! WD contains an array of diagonal matrices if (any(wd(:, 1, :) <= zero)) then info = max(30001, info + 1) return end if else ! WD contains an array of full M by M positive definite matrices do i = 1, n do j = 1, m wrk4(1:j, j) = wd(i, 1:j, j) end do call fctr(.false., wrk4, m, m, finfo) if (finfo /= 0) then info = max(30001, info + 1) return end if end do end if end if end if end subroutine fctrw