Scale t
by the inverse of scl
, i.e., compute t/scl
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n |
The number of rows of data in |
||
integer, | intent(in) | :: | m |
The number of columns of data in |
||
real(kind=wp), | intent(in) | :: | scl(ldscl,m) |
The scale values. |
||
integer, | intent(in) | :: | ldscl |
The leading dimension of array |
||
real(kind=wp), | intent(in) | :: | t(ldt,m) |
The array to be inversely scaled by |
||
integer, | intent(in) | :: | ldt |
The leading dimension of array |
||
real(kind=wp), | intent(out) | :: | sclt(ldsclt,m) |
The inversely scaled matrix. |
||
integer, | intent(in) | :: | ldsclt |
The leading dimension of array |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=wp), | public | :: | temp | ||||
integer, | public | :: | i | ||||
integer, | public | :: | j |
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