Changeset b4cf2b for src/LinearAlgebra
- Timestamp:
- Dec 4, 2010, 11:56:27 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:
- 8e9ce1
- Parents:
- bbf1bd
- git-author:
- Frederik Heber <heber@…> (11/16/10 09:45:32)
- git-committer:
- Frederik Heber <heber@…> (12/04/10 23:56:27)
- Location:
- src/LinearAlgebra
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
src/LinearAlgebra/MatrixContent.cpp
rbbf1bd rb4cf2b 30 30 #include <gsl/gsl_vector.h> 31 31 #include <cmath> 32 #include <cassert> 32 33 #include <iostream> 34 #include <set> 33 35 34 36 using namespace std; … … 201 203 const MatrixContent MatrixContent::operator*(const MatrixContent &rhs) const 202 204 { 203 gsl_matrix *res = gsl_matrix_alloc(rows, columns);205 gsl_matrix *res = gsl_matrix_alloc(rows, rhs.columns); 204 206 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, content, rhs.content, 0.0, res); 205 207 // gsl_matrix is taken over by constructor, hence no free … … 345 347 346 348 /** Transform the matrix to its eigenbasis ans return resulting eigenvalues. 349 * Note that we only return real-space part in case of non-symmetric matrix. 347 350 * \warn return vector has to be freed'd 351 * TODO: encapsulate return value in boost::shared_ptr or in VectorContent. 348 352 * \return gsl_vector pointer to vector of eigenvalues 349 353 */ 350 354 gsl_vector* MatrixContent::transformToEigenbasis() 351 355 { 352 ASSERT (rows == columns, 353 "MatrixContent::transformToEigenbasis() - only implemented for square matrices: "+toString(rows)+"!="+toString(columns)+"!"); 354 gsl_eigen_symmv_workspace *T = gsl_eigen_symmv_alloc(rows); 355 gsl_vector *eval = gsl_vector_alloc(rows); 356 gsl_matrix *evec = gsl_matrix_alloc(rows, rows); 357 gsl_eigen_symmv(content, eval, evec, T); 358 gsl_eigen_symmv_free(T); 359 gsl_matrix_memcpy(content, evec); 360 gsl_matrix_free(evec); 361 return eval; 356 if (rows == columns) { // symmetric 357 gsl_eigen_symmv_workspace *T = gsl_eigen_symmv_alloc(rows); 358 gsl_vector *eval = gsl_vector_alloc(rows); 359 gsl_matrix *evec = gsl_matrix_alloc(rows, rows); 360 gsl_eigen_symmv(content, eval, evec, T); 361 gsl_eigen_symmv_free(T); 362 gsl_matrix_memcpy(content, evec); 363 gsl_matrix_free(evec); 364 return eval; 365 } else { // non-symmetric 366 //ASSERT(false, "MatrixContent::transformToEigenbasis() - only implemented for square matrices: "+toString(rows)+"!="+toString(columns)+"!"); 367 368 // blow up gsl_matrix in content to square matrix, fill other components with zero 369 const size_t greaterDimension = rows > columns ? rows : columns; 370 gsl_matrix *content_square = gsl_matrix_alloc(greaterDimension, greaterDimension); 371 for (size_t i=0; i<greaterDimension; i++) { 372 for (size_t j=0; j<greaterDimension; j++) { 373 const double value = ((i < rows) && (j < columns)) ? gsl_matrix_get(content,i,j) : 0.; 374 gsl_matrix_set(content_square, i,j, value); 375 } 376 } 377 378 // solve eigenvalue problem 379 gsl_eigen_nonsymmv_workspace *T = gsl_eigen_nonsymmv_alloc(rows); 380 gsl_vector_complex *eval = gsl_vector_complex_alloc(greaterDimension); 381 gsl_matrix_complex *evec = gsl_matrix_complex_alloc(greaterDimension, greaterDimension); 382 gsl_eigen_nonsymmv(content_square, eval, evec, T); 383 gsl_eigen_nonsymmv_free(T); 384 385 // copy eigenvectors real-parts into content_square and display 386 std::cout << "Resulting eigenvector matrix is ["; 387 for (size_t i=0; i<greaterDimension; i++) { 388 for (size_t j=0; j<greaterDimension; j++) { 389 std::cout << "(" << GSL_REAL(gsl_matrix_complex_get(evec,i,j)) 390 << "," << GSL_IMAG(gsl_matrix_complex_get(evec,i,j)) << ")"; 391 if (j < greaterDimension-1) 392 std::cout << " "; 393 gsl_matrix_set(content_square, i,j, GSL_REAL(gsl_matrix_complex_get(evec,i,j))); 394 } 395 if (i < greaterDimension-1) 396 std::cout << "; "; 397 } 398 std::cout << "]" << std::endl; 399 400 // scan for zero rows and columns to drop 401 std::set<size_t> RowDropList; 402 std::set<size_t> ColumnDropList; 403 for (size_t i=0; i<greaterDimension; i++) { // only copy real space part 404 double Rsum = 0.; 405 double Csum = 0.; 406 for (size_t j=0; j<greaterDimension; j++) { 407 Rsum += fabs(GSL_REAL(gsl_matrix_complex_get(evec,i,j))); 408 Csum += fabs(GSL_REAL(gsl_matrix_complex_get(evec,j,i))); 409 } 410 Rsum /= (double)greaterDimension; 411 Csum /= (double)greaterDimension; 412 if (Rsum < MYEPSILON) 413 RowDropList.insert(i); 414 if (Csum < MYEPSILON) 415 ColumnDropList.insert(i); 416 } 417 418 // copy real-parts of complex eigenvalues and eigenvectors 419 gsl_vector *eval_real = gsl_vector_alloc(greaterDimension); 420 size_t I=0; 421 size_t J=0; 422 for (size_t i=0; i<greaterDimension; i++) { // only copy real space part 423 if (RowDropList.find(i) == RowDropList.end()) { 424 for (size_t j=0; j<greaterDimension; j++) { 425 if (ColumnDropList.find(j) == ColumnDropList.end()) { 426 if (fabs(GSL_IMAG(gsl_matrix_complex_get(evec,i,j))) > MYEPSILON) 427 std::cerr << "MatrixContent::transformToEigenbasis() - WARNING: eigenvectors are complex-valued!" << std::endl; 428 gsl_matrix_set(content, I,J, GSL_REAL(gsl_matrix_complex_get(evec,i,j))); 429 J++; 430 } 431 } 432 if (fabs(GSL_IMAG(gsl_vector_complex_get(eval,I))) > MYEPSILON) 433 std::cerr << "MatrixContent::transformToEigenbasis() - WARNING: eigenvectors are complex-valued!" << std::endl; 434 gsl_vector_set(eval_real, I, GSL_REAL(gsl_vector_complex_get(eval, i))); 435 I++; 436 } 437 } 438 gsl_matrix_complex_free(evec); 439 gsl_vector_complex_free(eval); 440 return eval_real; 441 } 362 442 } 363 443 … … 479 559 return result; 480 560 } 561 562 std::ostream & operator<<(std::ostream &ost, const MatrixContent &mat) 563 { 564 ost << "["; 565 for (size_t i=0;i<mat.rows;i++) { 566 for (size_t j=0;j<mat.columns;j++) { 567 ost << mat.at(i,j); 568 if (j != mat.columns-1) 569 ost << " "; 570 } 571 if (i != mat.rows-1) 572 ost << "; "; 573 } 574 ost << "]"; 575 return ost; 576 } -
src/LinearAlgebra/MatrixContent.hpp
rbbf1bd rb4cf2b 35 35 #include <gsl/gsl_matrix.h> 36 36 37 #include <iosfwd> 38 37 39 class Vector; 38 40 … … 40 42 { 41 43 friend Vector operator*(const MatrixContent &mat,const Vector &vec); 44 friend std::ostream & operator<<(std::ostream &ost, const MatrixContent &mat); 42 45 friend class RealSpaceMatrix; 43 46 friend class LinearSystemOfEquations; … … 104 107 const MatrixContent operator*(const MatrixContent &mat,const double factor); 105 108 Vector operator*(const MatrixContent &mat,const Vector &vec); 109 std::ostream & operator<<(std::ostream &ost, const MatrixContent &mat); 106 110 107 111 #endif /* MATRIXCONTENT_HPP_ */
Note:
See TracChangeset
for help on using the changeset viewer.