- Timestamp:
- Apr 6, 2011, 1:50:40 PM (14 years ago)
- Branches:
- 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
- Children:
- cd406d
- Parents:
- 0e4aa4
- git-author:
- Frederik Heber <heber@…> (03/30/11 20:45:38)
- git-committer:
- Frederik Heber <heber@…> (04/06/11 13:50:40)
- Location:
- src/LinearAlgebra
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
src/LinearAlgebra/MatrixContent.cpp
r0e4aa4 r6f68d6 44 44 MatrixContent::MatrixContent(size_t _rows, size_t _columns) : 45 45 rows(_rows), 46 columns(_columns) 46 columns(_columns), 47 free_content_on_exit(true) 47 48 { 48 49 content = gsl_matrix_calloc(rows, columns); … … 58 59 MatrixContent::MatrixContent(size_t _rows, size_t _columns, MatrixBaseCase) : 59 60 rows(_rows), 60 columns(_columns) 61 columns(_columns), 62 free_content_on_exit(true) 61 63 {} 62 64 … … 68 70 MatrixContent::MatrixContent(size_t _rows, size_t _columns, const double *src) : 69 71 rows(_rows), 70 columns(_columns) 72 columns(_columns), 73 free_content_on_exit(true) 71 74 { 72 75 content = gsl_matrix_calloc(rows, columns); … … 84 87 } 85 88 89 /** Constructor that parses square matrix from a stream. 90 * 91 * \note Matrix dimensions can be preparsed via 92 * MatrixContent::preparseMatrixDimensions() without harming the stream. 93 * 94 * \param &inputstream stream to parse from 95 */ 96 MatrixContent::MatrixContent(size_t _row, size_t _column, std::istream &inputstream) : 97 content(NULL), 98 rows(_row), 99 columns(_column), 100 free_content_on_exit(true) 101 { 102 // allocate matrix and set contents 103 content = gsl_matrix_alloc(rows, columns); 104 105 size_t row = 0; 106 do { 107 std::string line; 108 getline(inputstream, line); 109 //std::cout << line << std::endl; 110 std::stringstream parseline(line); 111 // skip comments 112 if ((parseline.peek() == '#') || (parseline.str().empty())) 113 continue; 114 // break on empty lines 115 if (parseline.peek() == '\n') 116 break; 117 118 // parse line with values 119 std::vector<double> line_of_values; 120 do { 121 double value; 122 parseline >> value >> ws; 123 line_of_values.push_back(value); 124 } while (parseline.good()); 125 126 // check number of columns parsed 127 ASSERT(line_of_values.size() == columns, 128 "MatrixContent::MatrixContent() - row " 129 +toString(row)+" has a different number of columns " 130 +toString(line_of_values.size())+" than others before " 131 +toString(columns)+"."); 132 133 for (size_t column = 0; column < columns; ++column) 134 set(row, column, line_of_values[column]); 135 ++row; 136 } while (inputstream.good()); 137 // check number of rows parsed 138 ASSERT(row == rows, 139 "MatrixContent::MatrixContent() - insufficent number of rows " 140 +toString(row)+" < "+toString(rows)+" read."); 141 } 142 143 86 144 /** Constructor for class MatrixContent. 87 * We embed the given gls_matrix pointer within this class and set it to NULL 88 * afterwards. 89 * \param *src source gsl_matrix vector to embed within this class 145 * 146 * \param *src source gsl_matrix vector to copy and put into this class 90 147 */ 91 148 MatrixContent::MatrixContent(gsl_matrix *&src) : 92 149 rows(src->size1), 93 columns(src->size2) 150 columns(src->size2), 151 free_content_on_exit(true) 94 152 { 95 153 content = gsl_matrix_alloc(src->size1, src->size2); … … 104 162 MatrixContent::MatrixContent(const MatrixContent &src) : 105 163 rows(src.rows), 106 columns(src.columns) 164 columns(src.columns), 165 free_content_on_exit(true) 107 166 { 108 167 content = gsl_matrix_alloc(src.rows, src.columns); … … 115 174 MatrixContent::MatrixContent(const MatrixContent *src) : 116 175 rows(src->rows), 117 columns(src->columns) 176 columns(src->columns), 177 free_content_on_exit(true) 118 178 { 119 179 ASSERT(src != NULL, "MatrixContent::MatrixContent - pointer to source matrix is NULL!"); … … 126 186 MatrixContent::~MatrixContent() 127 187 { 128 gsl_matrix_free(content); 188 if (free_content_on_exit) 189 gsl_matrix_free(content); 129 190 } 130 191 … … 607 668 }; 608 669 670 /** Preparses a matrix in an input stream for row and column count. 671 * 672 * \note This does not change the get position within the stream. 673 * 674 * @param inputstream to parse matrix from 675 * @return pair with (row, column) 676 */ 677 std::pair<size_t, size_t> MatrixContent::preparseMatrixDimensions(std::istream &inputstream) 678 { 679 const std::streampos initial_position(inputstream.tellg()); 680 // std::cout << "We are at position " 681 // +toString((int)inputstream.tellg())+" from start of stream." << std::endl; 682 size_t rows = 0; 683 size_t columns = 0; 684 do { 685 std::string line; 686 getline(inputstream, line); 687 std::stringstream parseline(line); 688 // skip comments 689 if ((parseline.peek() == '#') || (line.empty())) 690 continue; 691 // break on empty lines 692 if (parseline.peek() == '\n') 693 break; 694 695 // parse line with values 696 std::vector<double> line_of_values; 697 do { 698 double value; 699 parseline >> value >> ws; 700 line_of_values.push_back(value); 701 } while (parseline.good()); 702 703 // set matrixdimension if first line 704 if (columns == 0) { 705 columns = line_of_values.size(); 706 } else { // otherwise check 707 ASSERT(columns == line_of_values.size(), 708 "MatrixContent::MatrixContent() - row " 709 +toString(rows)+" has a different number of columns " 710 +toString(line_of_values.size())+" than this matrix has " 711 +toString(columns)+"."); 712 } 713 ++rows; 714 } while (inputstream.good()); 715 // clear end of file flags 716 inputstream.clear(); 717 718 // reset to initial position 719 inputstream.seekg(initial_position); 720 // std::cout << "We are again at position " 721 // +toString((int)inputstream.tellg())+" from start of stream." << std::endl; 722 723 return make_pair(rows, columns); 724 } 725 726 /** Writes matrix content to ostream in such a manner that it can be re-parsed 727 * by the code in the constructor. 728 * 729 * @param ost stream to write to 730 */ 731 void MatrixContent::write(std::ostream &ost) const 732 { 733 for (size_t row = 0; row < rows; ++row) { 734 for (size_t column = 0; column < columns; ++column) 735 ost << at(row, column) << "\t"; 736 ost << std::endl; 737 } 738 } 739 740 609 741 /* ============================= Operators =============================== */ 610 742 -
src/LinearAlgebra/MatrixContent.hpp
r0e4aa4 r6f68d6 71 71 MatrixContent(size_t _rows, size_t _columns, MatrixBaseCase); 72 72 MatrixContent(size_t rows, size_t columns, const double *src); 73 MatrixContent(size_t rows, size_t columns, std::istream &inputstream); 73 74 MatrixContent(gsl_matrix *&src); 74 75 MatrixContent(const MatrixContent &src); … … 124 125 VectorContent *getDiagonalVector() const; 125 126 127 static std::pair<size_t, size_t> preparseMatrixDimensions(std::istream &inputstream); 128 void write(std::ostream &ost) const; 129 126 130 protected: 127 131 double *Pointer(size_t m, size_t n); … … 133 137 134 138 private: 139 bool free_content_on_exit; 135 140 }; 136 141 -
src/LinearAlgebra/unittests/MatrixContentUnitTest.cpp
r0e4aa4 r6f68d6 266 266 CPPUNIT_ASSERT_EQUAL( false, m->IsPositiveDefinite() ); 267 267 }; 268 269 270 /** Unit test for reading from and writing matrix to stream 271 * 272 */ 273 void MatrixContentTest::ReadWriteTest() 274 { 275 // set up matrix 276 for (int i=0;i<4;i++) 277 for (int j=0;j<3;j++) 278 m->set(i,j, i*3+j ); 279 280 // write to stream 281 std::stringstream matrixstream; 282 m->write(matrixstream); 283 284 // parse in dimensions and check 285 std::pair<size_t, size_t> matrixdimensions = MatrixContent::preparseMatrixDimensions(matrixstream); 286 CPPUNIT_ASSERT_EQUAL( (size_t)4, matrixdimensions.first ); 287 CPPUNIT_ASSERT_EQUAL( (size_t)3, matrixdimensions.second ); 288 // parse in matrix and check 289 MatrixContent* n = new MatrixContent(4,3, matrixstream); 290 CPPUNIT_ASSERT_EQUAL( *m, *n ); 291 292 // free matrix 293 delete n; 294 } -
src/LinearAlgebra/unittests/MatrixContentUnitTest.hpp
r0e4aa4 r6f68d6 29 29 CPPUNIT_TEST (ExchangeTest ); 30 30 CPPUNIT_TEST (PropertiesTest ); 31 CPPUNIT_TEST (ReadWriteTest ); 31 32 CPPUNIT_TEST_SUITE_END(); 32 33 … … 40 41 void ExchangeTest(); 41 42 void PropertiesTest(); 43 void ReadWriteTest(); 42 44 43 45 private:
Note:
See TracChangeset
for help on using the changeset viewer.