28 #ifndef _math_scmat_matrix_h 29 #define _math_scmat_matrix_h 36 #include <math/scmat/abstract.h> 42 class SymmSCMatrixdouble;
43 class DiagSCMatrixdouble;
45 class SCMatrixBlockIter;
46 class SCMatrixRectBlock;
47 class SCMatrixLTriBlock;
48 class SCMatrixDiagBlock;
49 class SCVectorSimpleBlock;
52 class RefSymmSCMatrix;
94 void set_element(
int i,
double val)
const;
95 void accumulate_element(
int i,
double val)
const;
96 double get_element(
int)
const;
102 double maxabs()
const;
105 void normalize()
const;
106 void randomize()
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;
121 void print(std::ostream&out)
const;
122 void print(
const char*title=0,
123 std::ostream&out=
ExEnv::out0(),
int precision=10)
const;
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);
195 void accumulate_row(
const RefSCVector&,
int)
const;
196 void accumulate_column(
const RefSCVector&,
int)
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;
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,
228 double trace()
const;
243 double determ()
const;
299 void set_element(
int,
int,
double)
const;
300 void accumulate_element(
int,
int,
double)
const;
301 double get_element(
int,
int)
const;
303 RefSCMatrix get_subblock(
int br,
int er,
int bc,
int ec);
305 void assign_subblock(
const RefSCMatrix&,
int br,
int er,
int bc,
int ec);
307 void accumulate_subblock(
const RefSCMatrix&,
int,
int,
int,
int);
313 void accumulate_symmetric_outer_product(
const RefSCVector&)
const;
315 void accumulate_symmetric_product(
const RefSCMatrix&)
const;
316 void accumulate_symmetric_sum(
const RefSCMatrix&)
const;
319 SCMatrix::Transform = SCMatrix::NormalTransform)
const;
321 SCMatrix::Transform = SCMatrix::NormalTransform)
const;
325 void randomize()
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;
340 double trace()
const;
344 void print(std::ostream&)
const;
345 void print(
const char*title=0,
354 double determ()
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;
420 void scale(
double)
const;
421 void assign(
double)
const;
422 void assign(
const double*)
const;
423 void convert(
double*)
const;
434 double trace()
const;
435 void print(std::ostream&)
const;
436 void print(
const char*title=0,
442 double determ()
const;
526 #ifdef INLINE_FUNCTIONS 527 #include <math/scmat/matrix_i.h> void restore(StateIn &)
Restores the matrix from StateIn object. The vector must have been initialized already.
The RefSCVector class is a smart pointer to an SCVector specialization.
Definition: matrix.h:55
RefSymmSCMatrix symmetric_outer_product() const
The outer product of this with itself is a symmetric matrix.
Serializes objects that derive from SavableState.
Definition: stateout.h:61
SCVectordouble operator()(int) const
Return an l-value that can be used to assign or retrieve an element.
RefSCMatrix outer_product(const RefSCVector &v) const
Return the outer product between this and v.
The SCMatrix class is the abstract base class for general double valued n by m matrices.
Definition: abstract.h:195
The RefDiagSCMatrix class is a smart pointer to an DiagSCMatrix specialization.
Definition: matrix.h:380
A template class that maintains references counts.
Definition: ref.h:332
SCVector & operator*() const
Returns a C++ reference to the reference counted object.
Definition: ref.h:390
The RefSymmSCMatrix class is a smart pointer to an SCSymmSCMatrix specialization. ...
Definition: matrix.h:261
RefSCVector operator-(const RefSCVector &a) const
Subtract two vectors.
The SymmSCMatrix class is the abstract base class for diagonal double valued matrices.
Definition: abstract.h:503
RefSCVector()
Initializes the vector pointer to 0.
The RefSCDimension class is a smart pointer to an SCDimension specialization.
Definition: dim.h:156
Restores objects that derive from SavableState.
Definition: statein.h:70
The SCVector class is the abstract base class for double valued vectors.
Definition: abstract.h:97
RefSCVector & operator=(SCVector *v)
Make this refer to v.
static std::ostream & out0()
Return an ostream that writes from node 0.
The RefSCMatrix class is a smart pointer to an SCMatrix specialization.
Definition: matrix.h:135
The SymmSCMatrix class is the abstract base class for symmetric double valued matrices.
Definition: abstract.h:364
SCVectordouble operator[](int) const
Return an l-value that can be used to assign or retrieve an element.
RefSCVector operator+(const RefSCVector &a) const
Add two vectors.