ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/math/SquareMatrix.hpp
(Generate patch)

Comparing trunk/OOPSE-3.0/src/math/SquareMatrix.hpp (file contents):
Revision 1567 by tim, Wed Oct 13 23:53:40 2004 UTC vs.
Revision 1569 by tim, Thu Oct 14 23:28:09 2004 UTC

# Line 83 | Line 83 | namespace oopse {
83               SquareMatrix<Real, Dim> result;
84  
85               return result;
86 <        }
86 >        }        
87  
88        
89
88          /** Returns the determinant of this matrix. */
89          double determinant() const {
90              double det;
# Line 113 | Line 111 | namespace oopse {
111              return true;
112          }
113  
114 <        /** Tests if this matrix is orthogona. */            
114 >        /** Tests if this matrix is orthogonal. */            
115          bool isOrthogonal() {
116              SquareMatrix<Real, Dim> tmp;
117  
118              tmp = *this * transpose();
119  
120 <            return tmp.isUnitMatrix();
120 >            return tmp.isDiagonal();
121          }
122  
123          /** Tests if this matrix is diagonal. */
# Line 144 | Line 142 | namespace oopse {
142              return true;
143          }        
144  
145 +        void diagonalize() {
146 +            jacobi(m, eigenValues, ortMat);
147 +        }
148 +
149 +        /**
150 +         * Finds the eigenvalues and eigenvectors of a symmetric matrix
151 +         * @param eigenvals a reference to a vector3 where the
152 +         * eigenvalues will be stored. The eigenvalues are ordered so
153 +         * that eigenvals[0] <= eigenvals[1] <= eigenvals[2].
154 +         * @return an orthogonal matrix whose ith column is an
155 +         * eigenvector for the eigenvalue eigenvals[i]
156 +         */
157 +        SquareMatrix<Real, Dim>  findEigenvectors(Vector<Real, Dim>& eigenValues) {
158 +            SquareMatrix<Real, Dim> ortMat;
159 +            
160 +            if ( !isSymmetric()){
161 +                throw();
162 +            }
163 +            
164 +            SquareMatrix<Real, Dim> m(*this);
165 +            jacobi(m, eigenValues, ortMat);
166 +
167 +            return ortMat;
168 +        }
169 +        /**
170 +         * Jacobi iteration routines for computing eigenvalues/eigenvectors of
171 +         * real symmetric matrix
172 +         *
173 +         * @return true if success, otherwise return false
174 +         * @param a source matrix
175 +         * @param w output eigenvalues
176 +         * @param v output eigenvectors
177 +         */
178 +        void jacobi(const SquareMatrix<Real, Dim>& a,
179 +                              Vector<Real, Dim>& w,
180 +                              SquareMatrix<Real, Dim>& v);
181      };//end SquareMatrix
182  
183 +
184 + #define ROT(a,i,j,k,l) g=a(i, j);h=a(k, l);a(i, j)=g-s*(h+g*tau);a(k, l)=h+s*(g-h*tau)
185 + #define MAX_ROTATIONS 60
186 +
187 + template<Real, int Dim>
188 + void SquareMatrix<Real, int Dim>::jacobi(SquareMatrix<Real, Dim>& a,
189 +                                                                       Vector<Real, Dim>& w,
190 +                                                                       SquareMatrix<Real, Dim>& v) {
191 +    const int N = Dim;                                                                      
192 +    int i, j, k, iq, ip;
193 +    double tresh, theta, tau, t, sm, s, h, g, c;
194 +    double tmp;
195 +    Vector<Real, Dim> b, z;
196 +
197 +    // initialize
198 +    for (ip=0; ip<N; ip++)
199 +    {
200 +        for (iq=0; iq<N; iq++) v(ip, iq) = 0.0;
201 +        v(ip, ip) = 1.0;
202 +    }
203 +    for (ip=0; ip<N; ip++)
204 +    {
205 +        b(ip) = w(ip) = a(ip, ip);
206 +        z(ip) = 0.0;
207 +    }
208 +
209 +    // begin rotation sequence
210 +    for (i=0; i<MAX_ROTATIONS; i++)
211 +    {
212 +        sm = 0.0;
213 +        for (ip=0; ip<2; ip++)
214 +        {
215 +            for (iq=ip+1; iq<N; iq++) sm += fabs(a(ip, iq));
216 +        }
217 +        if (sm == 0.0) break;
218 +
219 +        if (i < 4) tresh = 0.2*sm/(9);
220 +        else tresh = 0.0;
221 +
222 +        for (ip=0; ip<2; ip++)
223 +        {
224 +            for (iq=ip+1; iq<N; iq++)
225 +            {
226 +                g = 100.0*fabs(a(ip, iq));
227 +                if (i > 4 && (fabs(w(ip))+g) == fabs(w(ip))
228 +                    && (fabs(w(iq))+g) == fabs(w(iq)))
229 +                {
230 +                    a(ip, iq) = 0.0;
231 +                }
232 +                else if (fabs(a(ip, iq)) > tresh)
233 +                {
234 +                    h = w(iq) - w(ip);
235 +                    if ( (fabs(h)+g) == fabs(h)) t = (a(ip, iq)) / h;
236 +                    else
237 +                    {
238 +                        theta = 0.5*h / (a(ip, iq));
239 +                        t = 1.0 / (fabs(theta)+sqrt(1.0+theta*theta));
240 +                        if (theta < 0.0) t = -t;
241 +                    }
242 +                    c = 1.0 / sqrt(1+t*t);
243 +                    s = t*c;
244 +                    tau = s/(1.0+c);
245 +                    h = t*a(ip, iq);
246 +                    z(ip) -= h;
247 +                    z(iq) += h;
248 +                    w(ip) -= h;
249 +                    w(iq) += h;
250 +                    a(ip, iq)=0.0;
251 +                    for (j=0;j<ip-1;j++)
252 +                    {
253 +                        ROT(a,j,ip,j,iq);
254 +                    }
255 +                    for (j=ip+1;j<iq-1;j++)
256 +                    {
257 +                        ROT(a,ip,j,j,iq);
258 +                    }
259 +                    for (j=iq+1; j<N; j++)
260 +                    {
261 +                        ROT(a,ip,j,iq,j);
262 +                    }
263 +                    for (j=0; j<N; j++)
264 +                    {
265 +                        ROT(v,j,ip,j,iq);
266 +                    }
267 +                }
268 +            }
269 +        }
270 +
271 +        for (ip=0; ip<N; ip++)
272 +        {
273 +            b(ip) += z(ip);
274 +            w(ip) = b(ip);
275 +            z(ip) = 0.0;
276 +        }
277 +    }
278 +
279 +    if ( i >= MAX_ROTATIONS )
280 +        return false;
281 +
282 +    // sort eigenfunctions
283 +    for (j=0; j<N; j++)
284 +    {
285 +        k = j;
286 +        tmp = w(k);
287 +        for (i=j; i<N; i++)
288 +        {
289 +            if (w(i) >= tmp)
290 +            {
291 +                k = i;
292 +                tmp = w(k);
293 +            }
294 +        }
295 +        if (k != j)
296 +        {
297 +            w(k) = w(j);
298 +            w(j) = tmp;
299 +            for (i=0; i<N; i++)
300 +            {
301 +                tmp = v(i, j);
302 +                v(i, j) = v(i, k);
303 +                v(i, k) = tmp;
304 +            }
305 +        }
306 +    }
307 +
308 +    //    insure eigenvector consistency (i.e., Jacobi can compute
309 +    //    vectors that are negative of one another (.707,.707,0) and
310 +    //    (-.707,-.707,0). This can reek havoc in
311 +    //    hyperstreamline/other stuff. We will select the most
312 +    //    positive eigenvector.
313 +    int numPos;
314 +    for (j=0; j<N; j++)
315 +    {
316 +        for (numPos=0, i=0; i<N; i++) if ( v(i, j) >= 0.0 ) numPos++;
317 +        if ( numPos < 2 ) for(i=0; i<N; i++) v(i, j) *= -1.0;
318 +    }
319 +
320 +    return true;
321   }
322 +
323 + #undef ROT
324 + #undef MAX_ROTATIONS
325 +
326 + }
327 +
328 +
329 + }
330   #endif //MATH_SQUAREMATRIX_HPP

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines