scale_vec Subroutine

public pure subroutine scale_vec(n, m, scl, ldscl, t, ldt, sclt, ldsclt)

Uses

  • proc~~scale_vec~~UsesGraph proc~scale_vec scale_vec module~odrpack_kinds odrpack_kinds proc~scale_vec->module~odrpack_kinds iso_fortran_env iso_fortran_env module~odrpack_kinds->iso_fortran_env

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.


Called by

proc~~scale_vec~~CalledByGraph proc~scale_vec scale_vec proc~trust_region trust_region proc~trust_region->proc~scale_vec proc~odr odr proc~odr->proc~trust_region 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 :: j

Source Code

   pure subroutine scale_vec(n, m, scl, ldscl, t, ldt, sclt, ldsclt)
   !! Scale `t` by the inverse of `scl`, i.e., compute `t/scl`.

      use odrpack_kinds, only: zero

      integer, intent(in) :: n
         !! Number of rows of data in `t`.
      integer, intent(in) :: m
         !! Number of columns of data in `t`.
      real(wp), intent(in) :: scl(ldscl, m)
         !! Scale values.
      integer, intent(in) :: ldscl
         !! Leading dimension of array `scl`.
      real(wp), intent(in) :: t(ldt, m)
         !! Array to be inversely scaled by `scl`.
      integer, intent(in) :: ldt
         !! Leading dimension of array `t`.
      real(wp), intent(out) :: sclt(ldsclt, m)
         !! Inversely scaled matrix.
      integer, intent(in) :: ldsclt
         !! Leading dimension of array `sclt`.

      ! Local scalars
      integer :: j

      ! Variable Definitions (alphabetically)
      !  J:       An indexing variable.

      if (n == 0 .or. m == 0) return

      if (scl(1, 1) >= zero) then
         if (ldscl >= n) then
            sclt(1:n, 1:m) = t(1:n, 1:m)/scl(1:n, 1:m)
         else
            do j = 1, m
               sclt(1:n, j) = t(1:n, j)/scl(1, j)
            end do
         end if
      else
         sclt(1:n, 1:m) = t(1:n, 1:m)/abs(scl(1, 1))
      end if

   end subroutine scale_vec