1 | //
|
---|
2 | // matrix.h
|
---|
3 | //
|
---|
4 | // Copyright (C) 1996 Limit Point Systems, Inc.
|
---|
5 | //
|
---|
6 | // Author: Curtis Janssen <cljanss@limitpt.com>
|
---|
7 | // Maintainer: LPS
|
---|
8 | //
|
---|
9 | // This file is part of the SC Toolkit.
|
---|
10 | //
|
---|
11 | // The SC Toolkit is free software; you can redistribute it and/or modify
|
---|
12 | // it under the terms of the GNU Library General Public License as published by
|
---|
13 | // the Free Software Foundation; either version 2, or (at your option)
|
---|
14 | // any later version.
|
---|
15 | //
|
---|
16 | // The SC Toolkit is distributed in the hope that it will be useful,
|
---|
17 | // but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
18 | // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
19 | // GNU Library General Public License for more details.
|
---|
20 | //
|
---|
21 | // You should have received a copy of the GNU Library General Public License
|
---|
22 | // along with the SC Toolkit; see the file COPYING.LIB. If not, write to
|
---|
23 | // the Free Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA.
|
---|
24 | //
|
---|
25 | // The U.S. Government is granted a limited license as per AL 91-7.
|
---|
26 | //
|
---|
27 |
|
---|
28 | #ifndef _math_scmat_matrix_h
|
---|
29 | #define _math_scmat_matrix_h
|
---|
30 | #ifdef __GNUC__
|
---|
31 | #pragma interface
|
---|
32 | #endif
|
---|
33 |
|
---|
34 | #include <iostream>
|
---|
35 |
|
---|
36 | #include <math/scmat/abstract.h>
|
---|
37 |
|
---|
38 | namespace sc {
|
---|
39 |
|
---|
40 | class SCVectordouble;
|
---|
41 | class SCMatrixdouble;
|
---|
42 | class SymmSCMatrixdouble;
|
---|
43 | class DiagSCMatrixdouble;
|
---|
44 |
|
---|
45 | class SCMatrixBlockIter;
|
---|
46 | class SCMatrixRectBlock;
|
---|
47 | class SCMatrixLTriBlock;
|
---|
48 | class SCMatrixDiagBlock;
|
---|
49 | class SCVectorSimpleBlock;
|
---|
50 |
|
---|
51 | class RefSCMatrix;
|
---|
52 | class RefSymmSCMatrix;
|
---|
53 | /** The RefSCVector class is a smart pointer to an SCVector specialization.
|
---|
54 | */
|
---|
55 | class RefSCVector: public Ref<SCVector> {
|
---|
56 | // standard overrides
|
---|
57 | public:
|
---|
58 | /** Initializes the vector pointer to 0. The reference must
|
---|
59 | be initialized before it is used. */
|
---|
60 | RefSCVector();
|
---|
61 | /// Make this and v refer to the same SCVector.
|
---|
62 | RefSCVector(const RefSCVector& v);
|
---|
63 | /// Make this refer to v.
|
---|
64 | RefSCVector(SCVector *v);
|
---|
65 | // don't allow automatic conversion from any reference to a
|
---|
66 | // described class
|
---|
67 | ~RefSCVector();
|
---|
68 | /// Make this refer to v.
|
---|
69 | RefSCVector& operator=(SCVector* v);
|
---|
70 | /// Make this and v refer to the same SCVector.
|
---|
71 | RefSCVector& operator=(const RefSCVector& v);
|
---|
72 |
|
---|
73 | // vector specific members
|
---|
74 | public:
|
---|
75 | /** Create a vector with dimension dim. The data values
|
---|
76 | are undefined. */
|
---|
77 | RefSCVector(const RefSCDimension& dim,const Ref<SCMatrixKit>&);
|
---|
78 |
|
---|
79 | /// Return an l-value that can be used to assign or retrieve an element.
|
---|
80 | SCVectordouble operator()(int) const;
|
---|
81 | /// Return an l-value that can be used to assign or retrieve an element.
|
---|
82 | SCVectordouble operator[](int) const;
|
---|
83 | /// Add two vectors.
|
---|
84 | RefSCVector operator+(const RefSCVector&a) const;
|
---|
85 | /// Subtract two vectors.
|
---|
86 | RefSCVector operator-(const RefSCVector&a) const;
|
---|
87 | /// Scale a vector.
|
---|
88 | RefSCVector operator*(double) const;
|
---|
89 | /// Return the outer product between this and v.
|
---|
90 | RefSCMatrix outer_product(const RefSCVector& v) const;
|
---|
91 | /// The outer product of this with itself is a symmetric matrix.
|
---|
92 | RefSymmSCMatrix symmetric_outer_product() const;
|
---|
93 |
|
---|
94 | void set_element(int i,double val) const;
|
---|
95 | void accumulate_element(int i,double val) const;
|
---|
96 | double get_element(int) const;
|
---|
97 | int n() const;
|
---|
98 | RefSCDimension dim() const;
|
---|
99 | Ref<SCMatrixKit> kit() const;
|
---|
100 | RefSCVector clone() const;
|
---|
101 | RefSCVector copy() const;
|
---|
102 | double maxabs() const;
|
---|
103 | double scalar_product(const RefSCVector&) const;
|
---|
104 | double dot(const RefSCVector&) const;
|
---|
105 | void normalize() const;
|
---|
106 | void randomize() const;
|
---|
107 | void assign(const RefSCVector& v) const;
|
---|
108 | void assign(double val) const;
|
---|
109 | void assign(const double* v) const;
|
---|
110 | void convert(double*) const;
|
---|
111 | void scale(double val) const;
|
---|
112 | void accumulate(const RefSCVector& v) const;
|
---|
113 | void accumulate_product(const RefSymmSCMatrix&, const RefSCVector&);
|
---|
114 | void accumulate_product(const RefSCMatrix&, const RefSCVector&);
|
---|
115 | void element_op(const Ref<SCElementOp>& op) const;
|
---|
116 | void element_op(const Ref<SCElementOp2>&,
|
---|
117 | const RefSCVector&) const;
|
---|
118 | void element_op(const Ref<SCElementOp3>&,
|
---|
119 | const RefSCVector&,
|
---|
120 | const RefSCVector&) const;
|
---|
121 | void print(std::ostream&out) const;
|
---|
122 | void print(const char*title=0,
|
---|
123 | std::ostream&out=ExEnv::out0(), int precision=10) const;
|
---|
124 | void save(StateOut&);
|
---|
125 | /// Restores the matrix from StateIn object. The vector must have been initialized already.
|
---|
126 | void restore(StateIn&);
|
---|
127 | };
|
---|
128 | RefSCVector operator*(double,const RefSCVector&);
|
---|
129 |
|
---|
130 | class RefSymmSCMatrix;
|
---|
131 | class RefDiagSCMatrix;
|
---|
132 | /** The RefSCMatrix class is a smart pointer to an SCMatrix
|
---|
133 | specialization.
|
---|
134 | */
|
---|
135 | class RefSCMatrix: public Ref<SCMatrix> {
|
---|
136 | // standard overrides
|
---|
137 | public:
|
---|
138 | /** Initializes the matrix pointer to 0. The reference must
|
---|
139 | be initialized before it is used. */
|
---|
140 | RefSCMatrix();
|
---|
141 | /// Make this and m refer to the same SCMatrix.
|
---|
142 | RefSCMatrix(const RefSCMatrix& m);
|
---|
143 | /// Make this refer to m.
|
---|
144 | RefSCMatrix(SCMatrix* m);
|
---|
145 | ~RefSCMatrix();
|
---|
146 | /// Make this refer to m.
|
---|
147 | RefSCMatrix& operator=(SCMatrix* m);
|
---|
148 | /// Make this and m refer to the same matrix.
|
---|
149 | RefSCMatrix& operator=(const RefSCMatrix& m);
|
---|
150 |
|
---|
151 | // matrix specific members
|
---|
152 | public:
|
---|
153 | /** Create a vector with dimension d1 by d2.
|
---|
154 | The data values are undefined. */
|
---|
155 | RefSCMatrix(const RefSCDimension& d1,const RefSCDimension& d2,
|
---|
156 | const Ref<SCMatrixKit>&);
|
---|
157 |
|
---|
158 | /// Multiply this by a vector and return a vector.
|
---|
159 | RefSCVector operator*(const RefSCVector&) const;
|
---|
160 |
|
---|
161 | /// Multiply this by a matrix and return a matrix.
|
---|
162 | RefSCMatrix operator*(const RefSCMatrix&) const;
|
---|
163 | RefSCMatrix operator*(const RefSymmSCMatrix&) const;
|
---|
164 | RefSCMatrix operator*(const RefDiagSCMatrix&) const;
|
---|
165 |
|
---|
166 | /// Multiply this by a scalar and return the result.
|
---|
167 | RefSCMatrix operator*(double) const;
|
---|
168 |
|
---|
169 | /// Matrix addition.
|
---|
170 | RefSCMatrix operator+(const RefSCMatrix&) const;
|
---|
171 | /// Matrix subtraction.
|
---|
172 | RefSCMatrix operator-(const RefSCMatrix&) const;
|
---|
173 |
|
---|
174 | /// Return the transpose of this.
|
---|
175 | RefSCMatrix t() const;
|
---|
176 | /// Return the inverse of this.
|
---|
177 | RefSCMatrix i() const;
|
---|
178 | /// Return the generalized inverse of this.
|
---|
179 | RefSCMatrix gi() const;
|
---|
180 |
|
---|
181 | /** These call the SCMatrix members of the same name
|
---|
182 | after checking for references to 0. */
|
---|
183 | RefSCMatrix clone() const;
|
---|
184 | RefSCMatrix copy() const;
|
---|
185 |
|
---|
186 | RefSCMatrix get_subblock(int br, int er, int bc, int ec);
|
---|
187 | void assign_subblock(const RefSCMatrix&, int br, int er, int bc, int ec,
|
---|
188 | int source_br = 0, int source_bc = 0);
|
---|
189 | void accumulate_subblock(const RefSCMatrix&, int, int, int, int,
|
---|
190 | int source_br = 0, int source_bc = 0);
|
---|
191 | RefSCVector get_row(int) const;
|
---|
192 | RefSCVector get_column(int) const;
|
---|
193 | void assign_row(const RefSCVector&, int) const;
|
---|
194 | void assign_column(const RefSCVector&, int) const;
|
---|
195 | void accumulate_row(const RefSCVector&, int) const;
|
---|
196 | void accumulate_column(const RefSCVector&, int) const;
|
---|
197 |
|
---|
198 | void accumulate_outer_product(const RefSCVector&,const RefSCVector&) const;
|
---|
199 | void accumulate_product(const RefSCMatrix&,const RefSCMatrix&) const;
|
---|
200 | void assign(const RefSCMatrix&) const;
|
---|
201 | void scale(double) const;
|
---|
202 | void randomize() const;
|
---|
203 | void assign(double) const;
|
---|
204 | void assign(const double*) const;
|
---|
205 | void assign(const double**) const;
|
---|
206 | void convert(double*) const;
|
---|
207 | void convert(double**) const;
|
---|
208 | void accumulate(const RefSCMatrix&) const;
|
---|
209 | void accumulate(const RefSymmSCMatrix&) const;
|
---|
210 | void accumulate(const RefDiagSCMatrix&) const;
|
---|
211 | void element_op(const Ref<SCElementOp>&) const;
|
---|
212 | void element_op(const Ref<SCElementOp2>&,
|
---|
213 | const RefSCMatrix&) const;
|
---|
214 | void element_op(const Ref<SCElementOp3>&,
|
---|
215 | const RefSCMatrix&,
|
---|
216 | const RefSCMatrix&) const;
|
---|
217 | int nrow() const;
|
---|
218 | int ncol() const;
|
---|
219 | RefSCDimension rowdim() const;
|
---|
220 | RefSCDimension coldim() const;
|
---|
221 | Ref<SCMatrixKit> kit() const;
|
---|
222 | void set_element(int,int,double) const;
|
---|
223 | void accumulate_element(int,int,double) const;
|
---|
224 | double get_element(int,int) const;
|
---|
225 | void print(std::ostream&) const;
|
---|
226 | void print(const char*title=0,
|
---|
227 | std::ostream&out=ExEnv::out0(), int =10) const;
|
---|
228 | double trace() const;
|
---|
229 | void save(StateOut&);
|
---|
230 | /// Restores the matrix from StateIn object. The matrix must have been initialized already.
|
---|
231 | void restore(StateIn&);
|
---|
232 |
|
---|
233 | /** Compute the singular value decomposition, this = U sigma V.t().
|
---|
234 | The dimension of sigma is the smallest dimension of this. U, V,
|
---|
235 | and sigma must already have the correct dimensions and are
|
---|
236 | overwritten. */
|
---|
237 | void svd(const RefSCMatrix &U,
|
---|
238 | const RefDiagSCMatrix &sigma,
|
---|
239 | const RefSCMatrix &V);
|
---|
240 | /** Solves this x = v. Overwrites v with x. */
|
---|
241 | double solve_lin(const RefSCVector& v) const;
|
---|
242 | /// Returns the determinant of the referenced matrix.
|
---|
243 | double determ() const;
|
---|
244 | /// Assign and examine matrix elements.
|
---|
245 | SCMatrixdouble operator()(int i,int j) const;
|
---|
246 |
|
---|
247 | /** If this matrix is blocked return the number of blocks.
|
---|
248 | * Otherwise return 1.
|
---|
249 | */
|
---|
250 | int nblock() const;
|
---|
251 | /** If this matrix is blocked return block i.
|
---|
252 | * Otherwise return this as block 0.
|
---|
253 | */
|
---|
254 | RefSCMatrix block(int i) const;
|
---|
255 | };
|
---|
256 | /// Allow multiplication with a scalar on the left.
|
---|
257 | RefSCMatrix operator*(double,const RefSCMatrix&);
|
---|
258 |
|
---|
259 | /** The RefSymmSCMatrix class is a smart pointer to an SCSymmSCMatrix
|
---|
260 | specialization. */
|
---|
261 | class RefSymmSCMatrix: public Ref<SymmSCMatrix> {
|
---|
262 | // standard overrides
|
---|
263 | public:
|
---|
264 | /** Initializes the matrix pointer to 0. The reference must
|
---|
265 | be initialized before it is used. */
|
---|
266 | RefSymmSCMatrix();
|
---|
267 | /// Make this and m refer to the same SCMatrix.
|
---|
268 | RefSymmSCMatrix(const RefSymmSCMatrix& m);
|
---|
269 | /// Make this refer to m.
|
---|
270 | RefSymmSCMatrix(SymmSCMatrix *m);
|
---|
271 | ~RefSymmSCMatrix();
|
---|
272 | /// Make this refer to m.
|
---|
273 | RefSymmSCMatrix& operator=(SymmSCMatrix* m);
|
---|
274 | /// Make this and m refer to the same matrix.
|
---|
275 | RefSymmSCMatrix& operator=(const RefSymmSCMatrix& m);
|
---|
276 |
|
---|
277 | // matrix specific members
|
---|
278 | public:
|
---|
279 | /** Create a vector with dimension d by d.
|
---|
280 | The data values are undefined. */
|
---|
281 | RefSymmSCMatrix(const RefSCDimension& d,const Ref<SCMatrixKit>&);
|
---|
282 | /// Multiply this by a matrix and return a matrix.
|
---|
283 | RefSCMatrix operator*(const RefSCMatrix&) const;
|
---|
284 | RefSCMatrix operator*(const RefSymmSCMatrix&) const;
|
---|
285 | /// Multiply this by a vector and return a vector.
|
---|
286 | RefSCVector operator*(const RefSCVector&a) const;
|
---|
287 | RefSymmSCMatrix operator*(double) const;
|
---|
288 | /// Matrix addition and subtraction.
|
---|
289 | RefSymmSCMatrix operator+(const RefSymmSCMatrix&) const;
|
---|
290 | RefSymmSCMatrix operator-(const RefSymmSCMatrix&) const;
|
---|
291 | /// Return the inverse of this.
|
---|
292 | RefSymmSCMatrix i() const;
|
---|
293 | /// Return the generalized inverse of this.
|
---|
294 | RefSymmSCMatrix gi() const;
|
---|
295 | /** These call the SCMatrix members of the same name after checking for
|
---|
296 | references to 0. */
|
---|
297 | RefSymmSCMatrix clone() const;
|
---|
298 | RefSymmSCMatrix copy() const;
|
---|
299 | void set_element(int,int,double) const;
|
---|
300 | void accumulate_element(int,int,double) const;
|
---|
301 | double get_element(int,int) const;
|
---|
302 |
|
---|
303 | RefSCMatrix get_subblock(int br, int er, int bc, int ec);
|
---|
304 | RefSymmSCMatrix get_subblock(int br, int er);
|
---|
305 | void assign_subblock(const RefSCMatrix&, int br, int er, int bc, int ec);
|
---|
306 | void assign_subblock(const RefSymmSCMatrix&, int br, int er);
|
---|
307 | void accumulate_subblock(const RefSCMatrix&, int, int, int, int);
|
---|
308 | void accumulate_subblock(const RefSymmSCMatrix&, int, int);
|
---|
309 | RefSCVector get_row(int);
|
---|
310 | void assign_row(const RefSCVector&, int);
|
---|
311 | void accumulate_row(const RefSCVector&, int);
|
---|
312 |
|
---|
313 | void accumulate_symmetric_outer_product(const RefSCVector&) const;
|
---|
314 | double scalar_product(const RefSCVector&) const;
|
---|
315 | void accumulate_symmetric_product(const RefSCMatrix&) const;
|
---|
316 | void accumulate_symmetric_sum(const RefSCMatrix&) const;
|
---|
317 | /// Add a * b * a.t() to this.
|
---|
318 | void accumulate_transform(const RefSCMatrix&a,const RefSymmSCMatrix&b,
|
---|
319 | SCMatrix::Transform = SCMatrix::NormalTransform) const;
|
---|
320 | void accumulate_transform(const RefSCMatrix&a,const RefDiagSCMatrix&b,
|
---|
321 | SCMatrix::Transform = SCMatrix::NormalTransform) const;
|
---|
322 | void accumulate_transform(const RefSymmSCMatrix&a,
|
---|
323 | const RefSymmSCMatrix&b) const;
|
---|
324 |
|
---|
325 | void randomize() const;
|
---|
326 | void assign(const RefSymmSCMatrix&) const;
|
---|
327 | void scale(double) const;
|
---|
328 | void assign(double) const;
|
---|
329 | void assign(const double*) const;
|
---|
330 | void assign(const double**) const;
|
---|
331 | void convert(double*) const;
|
---|
332 | void convert(double**) const;
|
---|
333 | void accumulate(const RefSymmSCMatrix&) const;
|
---|
334 | void element_op(const Ref<SCElementOp>&) const;
|
---|
335 | void element_op(const Ref<SCElementOp2>&,
|
---|
336 | const RefSymmSCMatrix&) const;
|
---|
337 | void element_op(const Ref<SCElementOp3>&,
|
---|
338 | const RefSymmSCMatrix&,
|
---|
339 | const RefSymmSCMatrix&) const;
|
---|
340 | double trace() const;
|
---|
341 | int n() const;
|
---|
342 | RefSCDimension dim() const;
|
---|
343 | Ref<SCMatrixKit> kit() const;
|
---|
344 | void print(std::ostream&) const;
|
---|
345 | void print(const char*title=0,
|
---|
346 | std::ostream&out=ExEnv::out0(), int =10) const;
|
---|
347 | void save(StateOut&);
|
---|
348 | /// Restores the matrix from StateIn object. The matrix must have been initialized already.
|
---|
349 | void restore(StateIn&);
|
---|
350 |
|
---|
351 | /** Solves this x = v. Overwrites v with x. */
|
---|
352 | double solve_lin(const RefSCVector&) const;
|
---|
353 | /// Returns the determinant of the referenced matrix.
|
---|
354 | double determ() const;
|
---|
355 | /// Returns the eigenvalues of the reference matrix.
|
---|
356 | RefDiagSCMatrix eigvals() const;
|
---|
357 | /// Returns the eigenvectors of the reference matrix.
|
---|
358 | RefSCMatrix eigvecs() const;
|
---|
359 | /** Sets eigvals to the eigenvalues and eigvecs
|
---|
360 | to the eigenvalues and eigenvectors of the referenced matrix.
|
---|
361 | The result satisfies eigvecs * eigvals * eigvecs.t() = (*this). */
|
---|
362 | void diagonalize(const RefDiagSCMatrix& eigvals,
|
---|
363 | const RefSCMatrix& eigvecs) const;
|
---|
364 | /// Assign and examine matrix elements.
|
---|
365 | SymmSCMatrixdouble operator()(int i,int j) const;
|
---|
366 | /** If this matrix is blocked return the number of blocks.
|
---|
367 | * Otherwise return 1.
|
---|
368 | */
|
---|
369 | int nblock() const;
|
---|
370 | /** If this matrix is blocked return block i.
|
---|
371 | * Otherwise return this as block 0.
|
---|
372 | */
|
---|
373 | RefSymmSCMatrix block(int i) const;
|
---|
374 | };
|
---|
375 | /// Allow multiplication with a scalar on the left.
|
---|
376 | RefSymmSCMatrix operator*(double,const RefSymmSCMatrix&);
|
---|
377 |
|
---|
378 | /** The RefDiagSCMatrix class is a smart pointer to an DiagSCMatrix
|
---|
379 | specialization. */
|
---|
380 | class RefDiagSCMatrix: public Ref<DiagSCMatrix> {
|
---|
381 | // standard overrides
|
---|
382 | public:
|
---|
383 | /** Initializes the matrix pointer to 0. The reference must
|
---|
384 | be initialized before it is used. */
|
---|
385 | RefDiagSCMatrix();
|
---|
386 | /// Make this and m refer to the same SCMatrix.
|
---|
387 | RefDiagSCMatrix(const RefDiagSCMatrix& m);
|
---|
388 | /// Make this refer to m.
|
---|
389 | RefDiagSCMatrix(DiagSCMatrix *m);
|
---|
390 | ~RefDiagSCMatrix();
|
---|
391 | /// Make this refer to m.
|
---|
392 | RefDiagSCMatrix& operator=(DiagSCMatrix* m);
|
---|
393 | /// Make this and m refer to the same matrix.
|
---|
394 | RefDiagSCMatrix& operator=(const RefDiagSCMatrix & m);
|
---|
395 |
|
---|
396 | // matrix specific members
|
---|
397 | public:
|
---|
398 | /** Create a diagonal matrix with dimension d by d. The data values
|
---|
399 | are undefined. */
|
---|
400 | RefDiagSCMatrix(const RefSCDimension&,const Ref<SCMatrixKit>&);
|
---|
401 | /// Multiply this by a matrix and return a matrix.
|
---|
402 | RefSCMatrix operator*(const RefSCMatrix&) const;
|
---|
403 | RefDiagSCMatrix operator*(double) const;
|
---|
404 | /// Matrix addition and subtraction.
|
---|
405 | RefDiagSCMatrix operator+(const RefDiagSCMatrix&) const;
|
---|
406 | RefDiagSCMatrix operator-(const RefDiagSCMatrix&) const;
|
---|
407 | /// Return the inverse of this.
|
---|
408 | RefDiagSCMatrix i() const;
|
---|
409 | /// Return the generalized inverse of this.
|
---|
410 | RefDiagSCMatrix gi() const;
|
---|
411 | /// These call the SCMatrix members of the same name
|
---|
412 | /// after checking for references to 0.
|
---|
413 | RefDiagSCMatrix clone() const;
|
---|
414 | RefDiagSCMatrix copy() const;
|
---|
415 | void set_element(int,double) const;
|
---|
416 | void accumulate_element(int,double) const;
|
---|
417 | double get_element(int) const;
|
---|
418 | void randomize() const;
|
---|
419 | void assign(const RefDiagSCMatrix&) const;
|
---|
420 | void scale(double) const;
|
---|
421 | void assign(double) const;
|
---|
422 | void assign(const double*) const;
|
---|
423 | void convert(double*) const;
|
---|
424 | void accumulate(const RefDiagSCMatrix&) const;
|
---|
425 | void element_op(const Ref<SCElementOp>&) const;
|
---|
426 | void element_op(const Ref<SCElementOp2>&,
|
---|
427 | const RefDiagSCMatrix&) const;
|
---|
428 | void element_op(const Ref<SCElementOp3>&,
|
---|
429 | const RefDiagSCMatrix&,
|
---|
430 | const RefDiagSCMatrix&) const;
|
---|
431 | int n() const;
|
---|
432 | RefSCDimension dim() const;
|
---|
433 | Ref<SCMatrixKit> kit() const;
|
---|
434 | double trace() const;
|
---|
435 | void print(std::ostream&) const;
|
---|
436 | void print(const char*title=0,
|
---|
437 | std::ostream&out=ExEnv::out0(), int =10) const;
|
---|
438 | void save(StateOut&);
|
---|
439 | /// Restores the matrix from StateIn object. The matrix must have been initialized already.
|
---|
440 | void restore(StateIn&);
|
---|
441 | /// Returns the determinant of the referenced matrix.
|
---|
442 | double determ() const;
|
---|
443 | /// Assign and examine matrix elements.
|
---|
444 | DiagSCMatrixdouble operator()(int i) const;
|
---|
445 | /** If this matrix is blocked return the number of blocks.
|
---|
446 | * Otherwise return 1.
|
---|
447 | */
|
---|
448 | int nblock() const;
|
---|
449 | /** If this matrix is blocked return block i.
|
---|
450 | * Otherwise return this as block 0.
|
---|
451 | */
|
---|
452 | RefDiagSCMatrix block(int i) const;
|
---|
453 | };
|
---|
454 | /// Allow multiplication with a scalar on the left.
|
---|
455 | RefDiagSCMatrix operator*(double,const RefDiagSCMatrix&);
|
---|
456 |
|
---|
457 | class SCVectordouble {
|
---|
458 | friend class RefSCVector;
|
---|
459 | private:
|
---|
460 | RefSCVector vector;
|
---|
461 | int i;
|
---|
462 |
|
---|
463 | SCVectordouble(SCVector*,int);
|
---|
464 | public:
|
---|
465 | SCVectordouble(const SCVectordouble&);
|
---|
466 | ~SCVectordouble();
|
---|
467 | double operator=(double a);
|
---|
468 | double operator=(const SCVectordouble&);
|
---|
469 | operator double();
|
---|
470 | double val() const;
|
---|
471 | };
|
---|
472 |
|
---|
473 | class SCMatrixdouble {
|
---|
474 | friend class RefSCMatrix;
|
---|
475 | private:
|
---|
476 | RefSCMatrix matrix;
|
---|
477 | int i;
|
---|
478 | int j;
|
---|
479 |
|
---|
480 | SCMatrixdouble(SCMatrix*,int,int);
|
---|
481 | public:
|
---|
482 | SCMatrixdouble(const SCMatrixdouble&);
|
---|
483 | ~SCMatrixdouble();
|
---|
484 | double operator=(double a);
|
---|
485 | double operator=(const SCMatrixdouble&);
|
---|
486 | operator double();
|
---|
487 | double val() const;
|
---|
488 | };
|
---|
489 |
|
---|
490 | class SymmSCMatrixdouble {
|
---|
491 | friend class RefSymmSCMatrix;
|
---|
492 | private:
|
---|
493 | RefSymmSCMatrix matrix;
|
---|
494 | int i;
|
---|
495 | int j;
|
---|
496 |
|
---|
497 | SymmSCMatrixdouble(SymmSCMatrix*,int,int);
|
---|
498 | public:
|
---|
499 | SymmSCMatrixdouble(const SCMatrixdouble&);
|
---|
500 | ~SymmSCMatrixdouble();
|
---|
501 | double operator=(double a);
|
---|
502 | double operator=(const SymmSCMatrixdouble&);
|
---|
503 | operator double();
|
---|
504 | double val() const;
|
---|
505 | };
|
---|
506 |
|
---|
507 | class DiagSCMatrixdouble {
|
---|
508 | friend class RefDiagSCMatrix;
|
---|
509 | private:
|
---|
510 | RefDiagSCMatrix matrix;
|
---|
511 | int i;
|
---|
512 | int j;
|
---|
513 |
|
---|
514 | DiagSCMatrixdouble(DiagSCMatrix*,int,int);
|
---|
515 | public:
|
---|
516 | DiagSCMatrixdouble(const SCMatrixdouble&);
|
---|
517 | ~DiagSCMatrixdouble();
|
---|
518 | double operator=(double a);
|
---|
519 | double operator=(const DiagSCMatrixdouble&);
|
---|
520 | operator double();
|
---|
521 | double val() const;
|
---|
522 | };
|
---|
523 |
|
---|
524 | }
|
---|
525 |
|
---|
526 | #ifdef INLINE_FUNCTIONS
|
---|
527 | #include <math/scmat/matrix_i.h>
|
---|
528 | #endif
|
---|
529 |
|
---|
530 | #endif
|
---|
531 |
|
---|
532 | // Local Variables:
|
---|
533 | // mode: c++
|
---|
534 | // c-file-style: "CLJ"
|
---|
535 | // End:
|
---|