DiagScale_spec.F
c**********************************************************************
#include "author.inc"
c* $Id: DiagScale_spec.F,v 1.1 1995/11/14 02:20:06 turner Exp $
c*
c* Apply diagonal (Jacobi) scaling.
c*
c* <PARAMETER LIST>
c*
c* Input:
c* idim - leading dimension of a
c* n - number of unknowns
c*
c* In/Out:
c* a - coefficient matrix in full format
c*
c* Output:
c* factor - vector of scaling factors
c* status - return status
c*
c* <SUBROUTINES REQUIRED>
c*
c* JT_FillVectorFloat
c* JT_y_eq_yx
c*
#include "copyright.inc"
c**********************************************************************
subroutine JT_DiagScale_Full (idim, n, a, factor, status)
implicit none
c
c ... Input:
integer idim, n
c
c ... In/Out:
real a(idim,n)
c
c ... Output:
integer status
real factor(n)
c
c ... Local:
integer i, j
real zero, one
c
parameter (zero=0.0d0, one=1.0d0)
c
c ... Compute and store inverse of diagonal elements.
call JT_FillVectorFloat (n, one, factor, status)
do i=1,n
if (a(i,i) .ne. zero) factor(i) = one / a(i,i)
enddo
c
c ... Scale coefficient.
do j=1,n
call JT_y_eq_yx (n, factor, a(1,j), status)
enddo
c
status = 0
return
end
c**********************************************************************
#include "author.inc"
c* $Id: DiagScale_spec.F,v 1.1 1995/11/14 02:20:06 turner Exp $
c*
c* Apply diagonal (Jacobi) scaling to a linear system with the
c* coefficient stored in ELL format.
c*
c* WARNING: This routine assumes that the diagonal is stored in
c* the first column.
c*
c* <PARAMETER LIST>
c*
c* Input:
c* idim - leading dimension of a
c* n - number of unknowns
c* maxnz - maximum number of non-zero elements in any row
c*
c* In/Out:
c* a - coefficient matrix in ELL format
c*
c* Output:
c* factor - vector of scaling factors
c* status - return status
c*
c* <SUBROUTINES REQUIRED>
c*
c* JT_FillVectorFloat
c* JT_y_eq_yx
c*
#include "copyright.inc"
c**********************************************************************
subroutine JT_DiagScale_ELL (idim, n, maxnz, a, factor, status)
implicit none
c
c ... Input:
integer idim, maxnz, n
c
c ... In/Out:
real a(idim,maxnz)
c
c ... Output:
integer status
real factor(n)
c
c ... Local:
integer i, j
real zero, one
c
parameter (zero=0.0d0, one=1.0d0)
c
c ... Compute and store inverse of diagonal elements.
call JT_FillVectorFloat (n, one, factor, status)
do i=1,n
if (a(i,1) .ne. zero) then
factor(i) = one / a(i,1)
a(i,1) = one
endif
enddo
c
c ... Scale coefficient.
do j=2,maxnz
call JT_y_eq_yx (n, factor, a(1,j), status)
enddo
c
status = 0
return
end
c**********************************************************************
#include "author.inc"
c* $Id: DiagScale_spec.F,v 1.1 1995/11/14 02:20:06 turner Exp $
c*
c* Apply diagonal (Jacobi) scaling to a linear system with the
c* coefficient stored in coordinate format.
c*
c* WARNING: This routine assumes that the diagonal is stored in
c* the first nelem elements.
c*
c* <PARAMETER LIST>
c*
c* Input:
c* nelem - number of active elements in a
c* n - number of unknowns
c* ja - index array
c*
c* In/Out:
c* a - coefficient matrix in COO format
c*
c* Output:
c* factor - vector of scaling factors
c* status - return status
c*
c* <SUBROUTINES REQUIRED>
c*
c* JT_FillVectorFloat
c* JT_y_eq_yx
c*
#include "copyright.inc"
c**********************************************************************
subroutine JT_DiagScale_COO (nelem, n, a, ja, factor, status)
implicit none
c
c ... Input:
integer nelem, n
integer ja(2*nelem)
c
c ... In/Out:
real a(nelem)
c
c ... Output:
integer status
real factor(n)
c
c ... Local:
integer i
real zero, one
c
parameter (zero=0.0d0, one=1.0d0)
c
c ... Compute and store inverse of diagonal elements.
call JT_FillVectorFloat (n, one, factor, status)
do i=1,n
if (a(i) .ne. zero) then
factor(i) = one / a(i)
a(i) = one
endif
enddo
c
c ... Scale coefficient.
do i=n+1,nelem
a(i) = factor(ja(i))*a(i)
enddo
c
status = 0
return
end
c**********************************************************************
#include "author.inc"
c* $Id: DiagScale_spec.F,v 1.1 1995/11/14 02:20:06 turner Exp $
c*
c* Apply diagonal (Jacobi) scaling to a linear system with the
c* coefficient stored in format.
c*
c* WARNING: This routine assumes that the diagonal is stored in
c* the first nelem elements.
c*
c* <PARAMETER LIST>
c*
c* Input:
c* nelem - number of active elements in a
c* n - number of unknowns
c* ja - index array
c*
c* In/Out:
c* a - coefficient matrix in RSS format
c*
c* Output:
c* factor - vector of scaling factors
c* status - return status
c*
c* <SUBROUTINES REQUIRED>
c*
c* JT_FillVectorFloat
c* JT_y_eq_yx
c*
#include "copyright.inc"
c**********************************************************************
subroutine JT_DiagScale_RSS (nelem, n, a, ja, factor, status)
implicit none
c
c ... Input:
integer nelem, n
integer ja(nelem)
c
c ... In/Out:
real a(nelem)
c
c ... Output:
integer status
real factor(n)
c
c ... Local:
integer i, ij
real zero, one
c
parameter (zero=0.0d0, one=1.0d0)
c
c ... Compute and store inverse of diagonal elements.
call JT_FillVectorFloat (n, one, factor, status)
do i=1,n
if (a(i) .ne. zero) then
factor(i) = one / a(i)
a(i) = one
endif
enddo
c
c ... Scale coefficient.
do ij=n+1,nelem
i = (ja(ij) - 1)/n + 1
a(ij) = factor(i)*a(ij)
enddo
c
status = 0
return
end
c**********************************************************************
#include "author.inc"
c* $Id: DiagScale_spec.F,v 1.1 1995/11/14 02:20:06 turner Exp $
c*
c* Apply diagonal (Jacobi) scaling to a linear system with the
c* coefficient stored in format.
c*
c* WARNING: This routine assumes that the diagonal is stored in
c* the first nelem elements.
c*
c* <PARAMETER LIST>
c*
c* Input:
c* nelem - number of active elements in a
c* n - number of unknowns
c* ja - index array
c*
c* In/Out:
c* a - coefficient matrix in CSS format
c*
c* Output:
c* factor - vector of scaling factors
c* status - return status
c*
c* <SUBROUTINES REQUIRED>
c*
c* JT_FillVectorFloat
c* JT_y_eq_yx
c*
#include "copyright.inc"
c**********************************************************************
subroutine JT_DiagScale_CSS (nelem, n, a, ja, factor, status)
implicit none
c
c ... Input:
integer nelem, n
integer ja(nelem)
c
c ... In/Out:
real a(nelem)
c
c ... Output:
integer status
real factor(n)
c
c ... Local:
integer i, j, ij
real zero, one
c
parameter (zero=0.0d0, one=1.0d0)
c
c ... Compute and store inverse of diagonal elements.
call JT_FillVectorFloat (n, one, factor, status)
do i=1,n
if (a(i) .ne. zero) then
factor(i) = one / a(i)
a(i) = one
endif
enddo
c
c ... Scale coefficient.
do ij=n+1,nelem
j = (ja(ij) - 1)/n + 1
i = ja(ij) - (j-1)*n
a(ij) = factor(i)*a(ij)
enddo
c
status = 0
return
end