loc_rwork Subroutine

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.


Called by

proc~~loc_rwork~~CalledByGraph proc~loc_rwork loc_rwork proc~access_workspace access_workspace proc~access_workspace->proc~loc_rwork proc~loc_rwork_c loc_rwork_c proc~loc_rwork_c->proc~loc_rwork proc~odr odr proc~odr->proc~loc_rwork proc~odr->proc~access_workspace proc~odr_long_c odr_long_c proc~odr_long_c->proc~odr proc~odr_medium_c odr_medium_c proc~odr_medium_c->proc~odr proc~odr_short_c odr_short_c proc~odr_short_c->proc~odr program~example1 example1 program~example1->proc~odr program~example2 example2 program~example2->proc~odr program~example3 example3 program~example3->proc~odr program~example4 example4 program~example4->proc~odr program~example5 example5 program~example5->proc~odr

Variables

Type Visibility Attributes Name Initial
integer, public :: next

Source Code

   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.

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

      ! Local scalars
      integer :: next

      ! Variable Definitions (alphabetically)
      !  NEXT:    The next available location with RWORK.

      if (n >= 1 .and. m >= 1 .and. np >= 1 .and. q >= 1 .and. ldwe >= 1 &
          .and. ld2we >= 1) then

         deltai = 1
         epsi = deltai + n*m
         xplusdi = epsi + n*q
         fni = xplusdi + n*m
         sdi = fni + n*q
         vcvi = sdi + np
         rvari = vcvi + np*np

         wssi = rvari + 1
         wssdeli = wssi + 1
         wssepsi = wssdeli + 1
         rcondi = wssepsi + 1
         etai = rcondi + 1
         olmavgi = etai + 1

         taui = olmavgi + 1
         alphai = taui + 1
         actrsi = alphai + 1
         pnormi = actrsi + 1
         rnormsi = pnormi + 1
         prersi = rnormsi + 1
         partoli = prersi + 1
         sstoli = partoli + 1
         taufaci = sstoli + 1
         epsmaci = taufaci + 1
         beta0i = epsmaci + 1

         betaci = beta0i + np
         betasi = betaci + np
         betani = betasi + np
         si = betani + np
         ssi = si + np
         ssfi = ssi + np
         qrauxi = ssfi + np
         ui = qrauxi + np
         fsi = ui + np

         fjacbi = fsi + n*q

         we1i = fjacbi + n*np*q

         diffi = we1i + ldwe*ld2we*q

         next = diffi + q*(np + m)

         if (isodr) then
            deltasi = next
            deltani = deltasi + n*m
            ti = deltani + n*m
            tti = ti + n*m
            omegai = tti + n*m
            fjacdi = omegai + q**2
            wrk1i = fjacdi + n*m*q
            next = wrk1i + n*m*q
         else
            deltasi = deltai
            deltani = deltai
            ti = deltai
            tti = deltai
            omegai = deltai
            fjacdi = deltai
            wrk1i = deltai
         end if

         wrk2i = next
         wrk3i = wrk2i + n*q
         wrk4i = wrk3i + np
         wrk5i = wrk4i + m*m
         wrk6i = wrk5i + m
         wrk7i = wrk6i + n*q*np
         loweri = wrk7i + 5*q
         upperi = loweri + np
         next = upperi + np

         lrwkmin = next
      else
         deltai = 1
         epsi = 1
         xplusdi = 1
         fni = 1
         sdi = 1
         vcvi = 1
         rvari = 1
         wssi = 1
         wssdeli = 1
         wssepsi = 1
         rcondi = 1
         etai = 1
         olmavgi = 1
         taui = 1
         alphai = 1
         actrsi = 1
         pnormi = 1
         rnormsi = 1
         prersi = 1
         partoli = 1
         sstoli = 1
         taufaci = 1
         epsmaci = 1
         beta0i = 1
         betaci = 1
         betasi = 1
         betani = 1
         si = 1
         ssi = 1
         ssfi = 1
         qrauxi = 1
         fsi = 1
         ui = 1
         fjacbi = 1
         we1i = 1
         diffi = 1
         deltasi = 1
         deltani = 1
         ti = 1
         tti = 1
         fjacdi = 1
         omegai = 1
         wrk1i = 1
         wrk2i = 1
         wrk3i = 1
         wrk4i = 1
         wrk5i = 1
         wrk6i = 1
         wrk7i = 1
         loweri = 1
         upperi = 1
         lrwkmin = 1
      end if

   end subroutine loc_rwork