jacobi.sa


Generated by gen_html_sa_files from ICSI. Contact gomes@icsi.berkeley.edu for details
 
------------------------->  GNU Sather - sourcefile  <-------------------------
-- Copyright (C) 2000 by K Hopper, University of Waikato, New Zealand        --
-- This file is part of the GNU Sather library. It is free software; you may --
-- redistribute  and/or modify it under the terms of the GNU Library General --
-- Public  License (LGPL)  as published  by the  Free  Software  Foundation; --
-- either version 2 of the license, or (at your option) any later version.   --
-- This  library  is distributed  in the  hope that it will  be  useful, but --
-- WITHOUT ANY WARRANTY without even the implied warranty of MERCHANTABILITY --
-- or FITNESS FOR A PARTICULAR PURPOSE. See Doc/LGPL for more details.       --
-- The license text is also available from:  Free Software Foundation, Inc., --
-- 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA                     --
-------------->  Please email comments to <bug-sather@gnu.org>  <--------------


class JACOBI

class JACOBI is -- This class implements the function which determines the eigenvalues -- and egienvectors of a real symmetric matrix. -- Version 1.2 March 2001. Copyright K Hopper, U of Waikato -- Development History -- ------------------- -- Date Who By Detail -- ---- ------ ------ -- 13 Jul 96 ac Original for Sather 1.1 Dist -- 5 Sep 97 kh Modified for cards/portability. -- 27 Mar 01 djw Completed and bugfixed. private const Max_Rotations : CARD := 99 ; private const Initial_Sweeps : CARD := 3 ; -- these use a non-zero threshold. private rotate( val : MAT, row1, col1, row2, col2 : CARD, inout elem1, inout elem2, scale, tau : FLT ) is -- This routine performs scaled rotations required when carrying out a -- Jacobi Transformation. elem1 := val[row1, col1] ; elem2 := val[row2, col2] ; val[row1, col1] := elem1 - scale * (elem2 + elem1 * tau) ; val[row2, col2] := elem2 + scale * (elem1 - elem2 * tau) end ; jacobi( src : MAT, out diagonal : VEC, out res : MAT ) : CARD pre src.size1 = src.size2 is -- This routine computes all eigenvalues and eigenvectors of the real -- symmetric matrix src. On output elements of src above the diagonal are -- destroyed. The diagonal argument returns the eigenvalues of src. The -- matrix res columns contain, on output, the normalized eigenvectors of src. -- The returned value is the number of Jacobi rotations which were required. coeff, elem1, elem2, scale, sum, tau, temp, thresh, theta : FLT ; loc_size : CARD := src.size1 ; loc_src : MAT := src.copy ; diagonal := VEC::create(loc_size) ; res := MAT::create(loc_size, loc_size) ; res.inplace(0.0); loop index : CARD := res.ind1! ; res[index, index] := 1.0 -- init to diag unit matrix end ; -- create and nitialize vectors vect1 and vect2 to have the -- values on the diagonal of loc_src. vect1 : VEC := VEC::create(loc_size) ; vect2 : VEC := VEC::create(loc_size) ; loop index : CARD := res.ind1! ; vect1[index] := loc_src[index, index] ; diagonal[index] := loc_src[index, index] ; vect2[index] := 0.0 end ; rot_cnt : CARD := 0 ; -- haven't done any rotns yet! loop outer : CARD := 0.upto!(Max_Rotations) ; -- count limit sum := 0.0 ; loop inner : CARD := 0.upto!(loc_size - 2) ; loop index : CARD := (inner + 1).upto!(loc_size - 1) ; sum := sum + loc_src[inner,index].abs end end ; -- Normal return, which relies on quadratic convergence to -- machine underflow if sum = 0.0 then -- ie no further change return rot_cnt end ; if outer < Initial_Sweeps then thresh := 0.2 * sum / (loc_size * loc_size).flt else -- the rest of the sweeps thresh := 0.0 end ; loop middle : CARD := 0.upto!(loc_size - 2) ; loop inner : CARD := (middle + 1).upto!(loc_size - 1) ; elem1 := 100.0 * loc_src[middle,inner].abs ; -- After 4 sweeps, skip rotation if the off diagonal -- element is small if (outer > Initial_Sweeps) and (diagonal[middle].abs + elem1 = diagonal[middle].abs) and (diagonal[inner].abs + elem1 = diagonal[inner].abs) then loc_src[middle,inner] := 0.0 elsif loc_src[middle,inner].abs > thresh then elem2 := diagonal[inner] - diagonal[middle] ; if elem2.abs + elem1 = elem2.abs then temp := loc_src[middle,inner] / elem2 else theta := 0.5 * elem2 / loc_src[middle,inner] ; temp := 1.0/(theta.abs + (1.0 + theta.square).sqrt) ; if theta < 0.0 then temp := -temp end end ; coeff := 1.0/(1.0 + temp * temp).sqrt ; scale := temp * coeff ; tau := scale / (1.0 + coeff) ; elem2 := temp * loc_src[middle,inner] ; vect2[middle] := vect2[middle] - elem2 ; vect2[inner] := vect2[inner] + elem2 ; diagonal[middle] := diagonal[middle] - elem2 ; diagonal[inner] := diagonal[inner] + elem2 ; loc_src[middle,inner] := 0.0 ; if middle >= 2 then loop index : CARD := 0.upto!(middle - 2) ; rotate(loc_src, index, middle, index, inner, inout elem2, inout elem1, scale, tau); end end ; if inner >= 2 then loop index : CARD := (middle + 1).upto!(inner - 2) ; rotate(loc_src, middle, index, index, inner, inout elem2, inout elem1, scale, tau) end end ; loop index : CARD := (inner + 1).upto!(loc_size - 1) ; rotate(loc_src, middle, index, inner, index, inout elem2, inout elem1, scale, tau) end ; loop index : CARD := 0.upto!(loc_size - 1) ; rotate(res, index, middle, index, inner, inout elem2, inout elem1, scale, tau) end ; rot_cnt := rot_cnt + 1 end end end ; loop inner : CARD := 0.upto!(loc_size - 1) ; vect1[inner] := vect1[inner] + vect2[inner] ; diagonal[inner] := vect1[inner] ; vect2[inner] := 0.0 end end ; raise "\nFatal error: too many iterations in jacobi\n" ; end ; -- jacobi end; -- class JACOBI