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

Comparing trunk/OOPSE-2.0/src/math/SquareMatrix.hpp (file contents):
Revision 1567 by tim, Wed Oct 13 23:53:40 2004 UTC vs.
Revision 1603 by tim, Tue Oct 19 21:28:55 2004 UTC

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

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines