dscale Subroutine

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

Uses

  • proc~~dscale~~UsesGraph proc~dscale dscale module~odrpack_kinds odrpack_kinds proc~dscale->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

The number of rows of data in t.

integer, intent(in) :: m

The number of columns of data in t.

real(kind=wp), intent(in) :: scl(ldscl,m)

The scale values.

integer, intent(in) :: ldscl

The leading dimension of array scl.

real(kind=wp), intent(in) :: t(ldt,m)

The array to be inversely scaled by scl.

integer, intent(in) :: ldt

The leading dimension of array t.

real(kind=wp), intent(out) :: sclt(ldsclt,m)

The inversely scaled matrix.

integer, intent(in) :: ldsclt

The leading dimension of array sclt.


Called by

proc~~dscale~~CalledByGraph proc~dscale dscale proc~dodlm dodlm proc~dodlm->proc~dscale proc~dodmn dodmn proc~dodmn->proc~dodlm proc~doddrv doddrv proc~doddrv->proc~dodmn proc~dodcnt dodcnt proc~dodcnt->proc~doddrv proc~odr odr proc~odr->proc~dodcnt 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
real(kind=wp), public :: temp
integer, public :: i
integer, public :: j

Source Code

   pure subroutine dscale(n, m, scl, ldscl, t, ldt, sclt, ldsclt)
   !! Scale `t` by the inverse of `scl`, i.e., compute `t/scl`.
      ! Routines Called (NONE)
      ! Date Written   860529   (YYMMDD)
      ! Revision Date  920304   (YYMMDD)

      use odrpack_kinds, only: zero, one

      integer, intent(in) :: n
         !! The number of rows of data in `t`.
      integer, intent(in) :: m
         !! The number of columns of data in `t`.
      real(wp), intent(in) :: scl(ldscl, m)
         !! The scale values.
      integer, intent(in) :: ldscl
         !! The leading dimension of array `scl`.
      real(wp), intent(in) :: t(ldt, m)
         !! The array to be inversely scaled by `scl`.
      integer, intent(in) :: ldt
         !! The leading dimension of array `t`.
      real(wp), intent(out) :: sclt(ldsclt, m)
         !! The inversely scaled matrix.
      integer, intent(in) :: ldsclt
         !! The leading dimension of array `sclt`.

      ! Local scalars
      real(wp) :: temp
      integer :: i, j

      ! Variable Definitions (alphabetically)
      !  I:       An indexing variable.
      !  J:       An indexing variable.
      !  LDSCL:   The leading dimension of array SCL.
      !  LDSCLT:  The leading dimension of array SCLT.
      !  LDT:     The leading dimension of array T.
      !  M:       The number of columns of data in T.
      !  N:       The number of rows of data in T.
      !  SCL:     The scale values.
      !  SCLT:    The inversely scaled matrix.
      !  T:       The array to be inversely scaled by SCL.
      !  TEMP:    A temporary scalar.

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

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

   end subroutine dscale