source: src/LinearAlgebra/MatrixContent.cpp@ 6d5a10

Action_Thermostats Add_AtomRandomPerturbation Add_FitFragmentPartialChargesAction Add_RotateAroundBondAction Add_SelectAtomByNameAction Added_ParseSaveFragmentResults AddingActions_SaveParseParticleParameters Adding_Graph_to_ChangeBondActions Adding_MD_integration_tests Adding_ParticleName_to_Atom Adding_StructOpt_integration_tests AtomFragments Automaking_mpqc_open AutomationFragmentation_failures Candidate_v1.5.4 Candidate_v1.6.0 Candidate_v1.6.1 ChangeBugEmailaddress ChangingTestPorts ChemicalSpaceEvaluator CombiningParticlePotentialParsing Combining_Subpackages Debian_Package_split Debian_package_split_molecuildergui_only Disabling_MemDebug Docu_Python_wait EmpiricalPotential_contain_HomologyGraph EmpiricalPotential_contain_HomologyGraph_documentation Enable_parallel_make_install Enhance_userguide Enhanced_StructuralOptimization Enhanced_StructuralOptimization_continued Example_ManyWaysToTranslateAtom Exclude_Hydrogens_annealWithBondGraph FitPartialCharges_GlobalError Fix_BoundInBox_CenterInBox_MoleculeActions Fix_ChargeSampling_PBC Fix_ChronosMutex Fix_FitPartialCharges Fix_FitPotential_needs_atomicnumbers Fix_ForceAnnealing Fix_IndependentFragmentGrids Fix_ParseParticles Fix_ParseParticles_split_forward_backward_Actions Fix_PopActions Fix_QtFragmentList_sorted_selection Fix_Restrictedkeyset_FragmentMolecule Fix_StatusMsg Fix_StepWorldTime_single_argument Fix_Verbose_Codepatterns Fix_fitting_potentials Fixes ForceAnnealing_goodresults ForceAnnealing_oldresults ForceAnnealing_tocheck ForceAnnealing_with_BondGraph ForceAnnealing_with_BondGraph_continued ForceAnnealing_with_BondGraph_continued_betteresults ForceAnnealing_with_BondGraph_contraction-expansion FragmentAction_writes_AtomFragments FragmentMolecule_checks_bonddegrees GeometryObjects Gui_Fixes Gui_displays_atomic_force_velocity ImplicitCharges IndependentFragmentGrids IndependentFragmentGrids_IndividualZeroInstances IndependentFragmentGrids_IntegrationTest IndependentFragmentGrids_Sole_NN_Calculation JobMarket_RobustOnKillsSegFaults JobMarket_StableWorkerPool JobMarket_unresolvable_hostname_fix MoreRobust_FragmentAutomation ODR_violation_mpqc_open PartialCharges_OrthogonalSummation PdbParser_setsAtomName PythonUI_with_named_parameters QtGui_reactivate_TimeChanged_changes Recreated_GuiChecks Rewrite_FitPartialCharges RotateToPrincipalAxisSystem_UndoRedo SaturateAtoms_findBestMatching SaturateAtoms_singleDegree StoppableMakroAction Subpackage_CodePatterns Subpackage_JobMarket Subpackage_LinearAlgebra Subpackage_levmar Subpackage_mpqc_open Subpackage_vmg Switchable_LogView ThirdParty_MPQC_rebuilt_buildsystem TrajectoryDependenant_MaxOrder TremoloParser_IncreasedPrecision TremoloParser_MultipleTimesteps TremoloParser_setsAtomName Ubuntu_1604_changes stable
Last change on this file since 6d5a10 was 6d5a10, checked in by Frederik Heber <heber@…>, 14 years ago

Added/Sorted some includes in LinearAlgebra and Exceptions.

  • Property mode set to 100644
File size: 18.1 KB
Line 
1/*
2 * MatrixContent.cpp
3 *
4 * Created on: Nov 14, 2010
5 * Author: heber
6 */
7
8
9// include config.h
10#ifdef HAVE_CONFIG_H
11#include <config.h>
12#endif
13
14#include "Helpers/MemDebug.hpp"
15
16#include "LinearAlgebra/RealSpaceMatrix.hpp"
17#include "Exceptions/NotInvertibleException.hpp"
18#include "Helpers/Assert.hpp"
19#include "Helpers/defs.hpp"
20#include "Helpers/fast_functions.hpp"
21#include "LinearAlgebra/Vector.hpp"
22#include "LinearAlgebra/VectorContent.hpp"
23#include "LinearAlgebra/MatrixContent.hpp"
24
25#include <gsl/gsl_blas.h>
26#include <gsl/gsl_eigen.h>
27#include <gsl/gsl_linalg.h>
28#include <gsl/gsl_matrix.h>
29#include <gsl/gsl_multimin.h>
30#include <gsl/gsl_vector.h>
31#include <cmath>
32#include <cassert>
33#include <iostream>
34#include <set>
35
36using namespace std;
37
38
39/** Constructor for class MatrixContent.
40 * \param rows number of rows
41 * \param columns number of columns
42 */
43MatrixContent::MatrixContent(size_t _rows, size_t _columns) :
44 rows(_rows),
45 columns(_columns)
46{
47 content = gsl_matrix_calloc(rows, columns);
48}
49
50/** Constructor of class VectorContent.
51 * We need this MatrixBaseCase for the VectorContentView class.
52 * There no content should be allocated, as it is just a view with an internal
53 * gsl_vector_view. Hence, MatrixBaseCase is just dummy class to give the
54 * constructor a unique signature.
55 * \param MatrixBaseCase
56 */
57MatrixContent::MatrixContent(size_t _rows, size_t _columns, MatrixBaseCase) :
58 rows(_rows),
59 columns(_columns)
60{}
61
62/** Constructor for class MatrixContent.
63 * \param rows number of rows
64 * \param columns number of columns
65 * \param *src array with components to initialize matrix with
66 */
67MatrixContent::MatrixContent(size_t _rows, size_t _columns, const double *src) :
68 rows(_rows),
69 columns(_columns)
70{
71 content = gsl_matrix_calloc(rows, columns);
72 set(0,0, src[0]);
73 set(1,0, src[1]);
74 set(2,0, src[2]);
75
76 set(0,1, src[3]);
77 set(1,1, src[4]);
78 set(2,1, src[5]);
79
80 set(0,2, src[6]);
81 set(1,2, src[7]);
82 set(2,2, src[8]);
83}
84
85/** Constructor for class MatrixContent.
86 * We embed the given gls_matrix pointer within this class and set it to NULL
87 * afterwards.
88 * \param *src source gsl_matrix vector to embed within this class
89 */
90MatrixContent::MatrixContent(gsl_matrix *&src) :
91 rows(src->size1),
92 columns(src->size2)
93{
94 content = gsl_matrix_alloc(src->size1, src->size2);
95 gsl_matrix_memcpy(content,src);
96// content = src;
97// src = NULL;
98}
99
100/** Copy constructor for class MatrixContent.
101 * \param &src reference to source MatrixContent
102 */
103MatrixContent::MatrixContent(const MatrixContent &src) :
104 rows(src.rows),
105 columns(src.columns)
106{
107 content = gsl_matrix_alloc(src.rows, src.columns);
108 gsl_matrix_memcpy(content,src.content);
109}
110
111/** Copy constructor for class MatrixContent.
112 * \param *src pointer to source MatrixContent
113 */
114MatrixContent::MatrixContent(const MatrixContent *src) :
115 rows(src->rows),
116 columns(src->columns)
117{
118 ASSERT(src != NULL, "MatrixContent::MatrixContent - pointer to source matrix is NULL!");
119 content = gsl_matrix_alloc(src->rows, src->columns);
120 gsl_matrix_memcpy(content,src->content);
121}
122
123/** Destructor for class MatrixContent.
124 */
125MatrixContent::~MatrixContent()
126{
127 gsl_matrix_free(content);
128}
129
130/** Set matrix to identity.
131 */
132void MatrixContent::setIdentity()
133{
134 for(int i=rows;i--;){
135 for(int j=columns;j--;){
136 set(i,j,i==j);
137 }
138 }
139}
140
141/** Set all matrix components to zero.
142 */
143void MatrixContent::setZero()
144{
145 for(int i=rows;i--;){
146 for(int j=columns;j--;){
147 set(i,j,0.);
148 }
149 }
150}
151
152/** Set all matrix components to a given value.
153 * \param _value value to set each component to
154 */
155void MatrixContent::setValue(double _value)
156{
157 for(int i=rows;i--;){
158 for(int j=columns;j--;){
159 set(i,j,_value);
160 }
161 }
162}
163
164/** Copy operator for MatrixContent with self-assignment check.
165 * \param &src matrix to compare to
166 * \return reference to this
167 */
168MatrixContent &MatrixContent::operator=(const MatrixContent &src)
169{
170 if(&src!=this){
171 gsl_matrix_memcpy(content,src.content);
172 }
173 return *this;
174}
175
176/** Addition operator.
177 * \param &rhs matrix to add
178 * \return reference to this
179 */
180const MatrixContent &MatrixContent::operator+=(const MatrixContent &rhs)
181{
182 gsl_matrix_add(content, rhs.content);
183 return *this;
184}
185
186/** Subtraction operator.
187 * \param &rhs matrix to subtract
188 * \return reference to this
189 */
190const MatrixContent &MatrixContent::operator-=(const MatrixContent &rhs)
191 {
192 gsl_matrix_sub(content, rhs.content);
193 return *this;
194}
195
196/** Multiplication operator.
197 * Note that here matrix have to have same dimensions.
198 * \param &rhs matrix to multiply with
199 * \return reference to this
200 */
201const MatrixContent &MatrixContent::operator*=(const MatrixContent &rhs)
202{
203 ASSERT(rows == rhs.rows,
204 "MatrixContent::operator*=() - row dimension differ: "+toString(rows)+" != "+toString(rhs.rows)+".");
205 ASSERT(columns == rhs.columns,
206 "MatrixContent::operator*=() - columns dimension differ: "+toString(columns)+" != "+toString(rhs.columns)+".");
207 (*this) = (*this)*rhs;
208 return *this;
209}
210
211/** Multiplication with copy operator.
212 * \param &rhs matrix to multiply with
213 * \return reference to newly allocated MatrixContent
214 */
215const MatrixContent MatrixContent::operator*(const MatrixContent &rhs) const
216{
217 gsl_matrix *res = gsl_matrix_alloc(rows, rhs.columns);
218 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, content, rhs.content, 0.0, res);
219 // gsl_matrix is taken over by constructor, hence no free
220 MatrixContent tmp(res);
221 gsl_matrix_free(res);
222 return tmp;
223}
224
225/* ========================== Accessing =============================== */
226
227/** Accessor for manipulating component (i,j).
228 * \param i row number
229 * \param j column number
230 * \return reference to component (i,j)
231 */
232double &MatrixContent::at(size_t i, size_t j)
233{
234 ASSERT((i>=0) && (i<rows),
235 "MatrixContent::at() - Index i="+toString(i)+" for Matrix access out of range [0,"+toString(rows)+"]");
236 ASSERT((j>=0) && (j<columns),
237 "MatrixContent::at() - Index j="+toString(j)+" for Matrix access out of range [0,"+toString(columns)+"]");
238 return *gsl_matrix_ptr (content, i, j);
239}
240
241/** Constant accessor for (value of) component (i,j).
242 * \param i row number
243 * \param j column number
244 * \return const component (i,j)
245 */
246const double MatrixContent::at(size_t i, size_t j) const
247{
248 ASSERT((i>=0) && (i<rows),
249 "MatrixContent::at() - Index i="+toString(i)+" for Matrix access out of range [0,"+toString(rows)+"]");
250 ASSERT((j>=0) && (j<columns),
251 "MatrixContent::at() - Index j="+toString(j)+" for Matrix access out of range [0,"+toString(columns)+"]");
252 return gsl_matrix_get(content, i, j);
253}
254
255/** These functions return a pointer to the \a m-th element of a matrix.
256 * If \a m or \a n lies outside the allowed range of 0 to MatrixContent::dimension-1 then the error handler is invoked and a null pointer is returned.
257 * \param m index
258 * \return pointer to \a m-th element
259 */
260double *MatrixContent::Pointer(size_t m, size_t n)
261{
262 return gsl_matrix_ptr (content, m, n);
263};
264
265/** These functions return a constant pointer to the \a m-th element of a matrix.
266 * If \a m or \a n lies outside the allowed range of 0 to MatrixContent::dimension-1 then the error handler is invoked and a null pointer is returned.
267 * \param m index
268 * \return const pointer to \a m-th element
269 */
270const double *MatrixContent::const_Pointer(size_t m, size_t n) const
271{
272 return gsl_matrix_const_ptr (content, m, n);
273};
274
275/* ========================== Initializing =============================== */
276
277/** Setter for component (i,j).
278 * \param i row numbr
279 * \param j column numnber
280 * \param value value to set componnt (i,j) to
281 */
282void MatrixContent::set(size_t i, size_t j, const double value)
283{
284 ASSERT((i>=0) && (i<rows),
285 "MatrixContent::set() - Index i="+toString(i)+" for Matrix access out of range [0,"+toString(rows)+"]");
286 ASSERT((j>=0) && (j<columns),
287 "MatrixContent::set() - Index j="+toString(j)+" for Matrix access out of range [0,"+toString(columns)+"]");
288 gsl_matrix_set(content,i,j,value);
289}
290
291/** This function sets the matrix from a double array.
292 * Creates a matrix view of the array and performs a memcopy.
293 * \param *x array of values (no dimension check is performed)
294 */
295void MatrixContent::setFromDoubleArray(double * x)
296{
297 gsl_matrix_view m = gsl_matrix_view_array (x, rows, columns);
298 gsl_matrix_memcpy (content, &m.matrix);
299};
300
301/* ====================== Exchanging elements ============================ */
302/** This function exchanges the \a i-th and \a j-th row of the matrix in-place.
303 * \param i i-th row to swap with ...
304 * \param j ... j-th row to swap against
305 */
306bool MatrixContent::SwapRows(size_t i, size_t j)
307{
308 return (gsl_matrix_swap_rows (content, i, j) == GSL_SUCCESS);
309};
310
311/** This function exchanges the \a i-th and \a j-th column of the matrix in-place.
312 * \param i i-th column to swap with ...
313 * \param j ... j-th column to swap against
314 */
315bool MatrixContent::SwapColumns(size_t i, size_t j)
316{
317 return (gsl_matrix_swap_columns (content, i, j) == GSL_SUCCESS);
318};
319
320/** This function exchanges the \a i-th row and \a j-th column of the matrix in-place.
321 * The matrix must be square for this operation to be possible.
322 * \param i i-th row to swap with ...
323 * \param j ... j-th column to swap against
324 */
325bool MatrixContent::SwapRowColumn(size_t i, size_t j)
326{
327 assert (rows == columns && "The matrix must be square for swapping row against column to be possible.");
328 return (gsl_matrix_swap_rowcol (content, i, j) == GSL_SUCCESS);
329};
330
331/** Return transposed matrix.
332 * \return new matrix that is transposed of this.
333 */
334MatrixContent MatrixContent::transpose() const
335{
336 gsl_matrix *res = gsl_matrix_alloc(columns, rows); // column and row dimensions exchanged!
337 gsl_matrix_transpose_memcpy(res, content);
338 MatrixContent newContent(res);
339 gsl_matrix_free(res);
340 return newContent;
341}
342
343/** Turn this matrix into its transposed.
344 * Note that this is only possible if rows == columns.
345 */
346MatrixContent &MatrixContent::transpose()
347{
348 ASSERT( rows == columns,
349 "MatrixContent::transpose() - cannot transpose onto itself as matrix not square: "+toString(rows)+"!="+toString(columns)+"!");
350 double tmp;
351 for (size_t i=0;i<rows;i++)
352 for (size_t j=i+1;j<rows;j++) {
353 tmp = at(j,i);
354 at(j,i) = at(i,j);
355 at(i,j) = tmp;
356 }
357 return *this;
358}
359
360/** Transform the matrix to its eigenbasis ans return resulting eigenvalues.
361 * Note that we only return real-space part in case of non-symmetric matrix.
362 * \warn return vector has to be freed'd
363 * TODO: encapsulate return value in boost::shared_ptr or in VectorContent.
364 * \return gsl_vector pointer to vector of eigenvalues
365 */
366gsl_vector* MatrixContent::transformToEigenbasis()
367{
368 if (rows == columns) { // symmetric
369 gsl_eigen_symmv_workspace *T = gsl_eigen_symmv_alloc(rows);
370 gsl_vector *eval = gsl_vector_alloc(rows);
371 gsl_matrix *evec = gsl_matrix_alloc(rows, rows);
372 gsl_eigen_symmv(content, eval, evec, T);
373 gsl_eigen_symmv_free(T);
374 gsl_matrix_memcpy(content, evec);
375 gsl_matrix_free(evec);
376 return eval;
377 } else { // non-symmetric
378 //ASSERT(false, "MatrixContent::transformToEigenbasis() - only implemented for square matrices: "+toString(rows)+"!="+toString(columns)+"!");
379
380 // blow up gsl_matrix in content to square matrix, fill other components with zero
381 const size_t greaterDimension = rows > columns ? rows : columns;
382 gsl_matrix *content_square = gsl_matrix_alloc(greaterDimension, greaterDimension);
383 for (size_t i=0; i<greaterDimension; i++) {
384 for (size_t j=0; j<greaterDimension; j++) {
385 const double value = ((i < rows) && (j < columns)) ? gsl_matrix_get(content,i,j) : 0.;
386 gsl_matrix_set(content_square, i,j, value);
387 }
388 }
389
390 // show squared matrix
391 MatrixContent *ContentSquare = new MatrixViewContent(gsl_matrix_submatrix(content,0,0,content->size1, content->size2));
392 std::cout << "The squared matrix is " << *ContentSquare << std::endl;
393 delete ContentSquare;
394
395 // solve eigenvalue problem
396 gsl_eigen_nonsymmv_workspace *T = gsl_eigen_nonsymmv_alloc(rows);
397 gsl_vector_complex *eval = gsl_vector_complex_alloc(greaterDimension);
398 gsl_matrix_complex *evec = gsl_matrix_complex_alloc(greaterDimension, greaterDimension);
399 gsl_eigen_nonsymmv(content_square, eval, evec, T);
400 gsl_eigen_nonsymmv_free(T);
401
402 // copy eigenvectors real-parts into content_square and ...
403 // ... show complex-valued eigenvector matrix
404 std::cout << "Resulting eigenvector matrix is [";
405 for (size_t i=0; i<greaterDimension; i++) {
406 for (size_t j=0; j<greaterDimension; j++) {
407 std::cout << "(" << GSL_REAL(gsl_matrix_complex_get(evec,i,j))
408 << "," << GSL_IMAG(gsl_matrix_complex_get(evec,i,j)) << ")";
409 if (j < greaterDimension-1)
410 std::cout << " ";
411 gsl_matrix_set(content_square, i,j, GSL_REAL(gsl_matrix_complex_get(evec,i,j)));
412 }
413 if (i < greaterDimension-1)
414 std::cout << "; ";
415 }
416 std::cout << "]" << std::endl;
417
418 // scan for zero rows and columns to drop
419 std::set<size_t> RowDropList;
420 std::set<size_t> ColumnDropList;
421 for (size_t i=0; i<greaterDimension; i++) { // only copy real space part
422 double Rsum = 0.;
423 double Csum = 0.;
424 for (size_t j=0; j<greaterDimension; j++) {
425 Rsum += fabs(GSL_REAL(gsl_matrix_complex_get(evec,i,j)));
426 Csum += fabs(GSL_REAL(gsl_matrix_complex_get(evec,j,i)));
427 }
428 Rsum /= (double)greaterDimension;
429 Csum /= (double)greaterDimension;
430 if (Rsum < MYEPSILON)
431 RowDropList.insert(i);
432 if (Csum < MYEPSILON)
433 ColumnDropList.insert(i);
434 }
435
436 // copy real-parts of complex eigenvalues and eigenvectors
437 gsl_vector *eval_real = gsl_vector_alloc(greaterDimension);
438 size_t I=0;
439 size_t J=0;
440 for (size_t i=0; i<greaterDimension; i++) { // only copy real space part
441 if (RowDropList.find(i) == RowDropList.end()) {
442 for (size_t j=0; j<greaterDimension; j++) {
443 if (ColumnDropList.find(j) == ColumnDropList.end()) {
444 if (fabs(GSL_IMAG(gsl_matrix_complex_get(evec,i,j))) > MYEPSILON)
445 std::cerr << "MatrixContent::transformToEigenbasis() - WARNING: eigenvectors are complex-valued!" << std::endl;
446 gsl_matrix_set(content, I,J, GSL_REAL(gsl_matrix_complex_get(evec,i,j)));
447 J++;
448 }
449 }
450 if (fabs(GSL_IMAG(gsl_vector_complex_get(eval,I))) > MYEPSILON)
451 std::cerr << "MatrixContent::transformToEigenbasis() - WARNING: eigenvectors are complex-valued!" << std::endl;
452 gsl_vector_set(eval_real, I, GSL_REAL(gsl_vector_complex_get(eval, i)));
453 I++;
454 }
455 }
456 gsl_matrix_complex_free(evec);
457 gsl_vector_complex_free(eval);
458 return eval_real;
459 }
460}
461
462/* ============================ Properties ============================== */
463/** Checks whether matrix' elements are strictly null.
464 * \return true - is null, false - else
465 */
466bool MatrixContent::IsNull() const
467{
468 return gsl_matrix_isnull (content);
469};
470
471/** Checks whether matrix' elements are strictly positive.
472 * \return true - is positive, false - else
473 */
474bool MatrixContent::IsPositive() const
475{
476 return gsl_matrix_ispos (content);
477};
478
479/** Checks whether matrix' elements are strictly negative.
480 * \return true - is negative, false - else
481 */
482bool MatrixContent::IsNegative() const
483{
484 return gsl_matrix_isneg (content);
485};
486
487/** Checks whether matrix' elements are strictly non-negative.
488 * \return true - is non-negative, false - else
489 */
490bool MatrixContent::IsNonNegative() const
491{
492 return gsl_matrix_isnonneg (content);
493};
494
495/** This function performs a Cholesky decomposition to determine whether matrix is positive definite.
496 * We check whether GSL returns GSL_EDOM as error, indicating that decomposition failed due to matrix not being positive-definite.
497 * \return true - matrix is positive-definite, false - else
498 */
499bool MatrixContent::IsPositiveDefinite() const
500{
501 if (rows != columns) // only possible for square matrices.
502 return false;
503 else
504 return (gsl_linalg_cholesky_decomp (content) != GSL_EDOM);
505};
506
507
508/** Calculates the determinant of the matrix.
509 * if matrix is square, uses LU decomposition.
510 */
511double MatrixContent::Determinant() const
512{
513 int signum = 0;
514 assert (rows == columns && "Determinant can only be calculated for square matrices.");
515 gsl_permutation *p = gsl_permutation_alloc(rows);
516 gsl_linalg_LU_decomp(content, p, &signum);
517 gsl_permutation_free(p);
518 return gsl_linalg_LU_det(content, signum);
519};
520
521/* ============================= Operators =============================== */
522
523/** Scalar multiplication operator.
524 * \param factor factor to scale with
525 */
526const MatrixContent &MatrixContent::operator*=(const double factor)
527{
528 gsl_matrix_scale(content, factor);
529 return *this;
530}
531
532/** Scalar multiplication and copy operator.
533 * \param factor factor to scale with
534 * \param &mat MatrixContent to scale
535 * \return copied and scaled MatrixContent
536 */
537const MatrixContent operator*(const double factor,const MatrixContent& mat)
538{
539 MatrixContent tmp = mat;
540 tmp*=factor;
541 return tmp;
542}
543
544/** Scalar multiplication and copy operator (with operands exchanged).
545 * \param &mat MatrixContent to scale
546 * \param factor factor to scale with
547 * \return copied and scaled MatrixContent
548 */
549const MatrixContent operator*(const MatrixContent &mat,const double factor)
550{
551 return factor*mat;
552}
553
554/** Equality operator.
555 * Note that we use numerical sensible checking, i.e. with threshold MYEPSILON.
556 * \param &rhs MatrixContent to checks against
557 */
558bool MatrixContent::operator==(const MatrixContent &rhs) const
559 {
560 if ((rows == rhs.rows) && (columns == rhs.columns)) {
561 for(int i=rows;i--;){
562 for(int j=columns;j--;){
563 if(fabs(at(i,j)-rhs.at(i,j))>MYEPSILON){
564 return false;
565 }
566 }
567 }
568 return true;
569 }
570 return false;
571}
572
573Vector operator*(const MatrixContent &mat,const Vector &vec)
574{
575 Vector result;
576 gsl_blas_dgemv( CblasNoTrans, 1.0, mat.content, vec.content->content, 0.0, result.content->content);
577 return result;
578}
579
580std::ostream & operator<<(std::ostream &ost, const MatrixContent &mat)
581{
582 ost << "[";
583 for (size_t i=0;i<mat.rows;i++) {
584 for (size_t j=0;j<mat.columns;j++) {
585 ost << mat.at(i,j);
586 if (j != mat.columns-1)
587 ost << " ";
588 }
589 if (i != mat.rows-1)
590 ost << "; ";
591 }
592 ost << "]";
593 return ost;
594}
Note: See TracBrowser for help on using the repository browser.