source: src/LinearAlgebra/gslmatrix.cpp@ bf3817

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 bf3817 was bf3817, checked in by Frederik Heber <heber@…>, 15 years ago

Added ifdef HAVE_CONFIG and config.h include to each and every cpp file.

  • is now topmost in front of MemDebug.hpp (and any other).
  • Property mode set to 100644
File size: 7.4 KB
Line 
1/*
2 * gslmatrix.cpp
3 *
4 * Created on: Jan 8, 2010
5 * Author: heber
6 */
7
8// include config.h
9#ifdef HAVE_CONFIG_H
10#include <config.h>
11#endif
12
13#include "Helpers/MemDebug.hpp"
14
15using namespace std;
16
17#include "LinearAlgebra/gslmatrix.hpp"
18#include "Helpers/helpers.hpp"
19#include "Helpers/fast_functions.hpp"
20
21#include <cassert>
22#include <gsl/gsl_linalg.h>
23
24/** Constructor of class GSLMatrix.
25 * Allocates GSL structures
26 * \param m dimension of matrix
27 */
28GSLMatrix::GSLMatrix(size_t m, size_t n) : rows(m), columns(n)
29{
30 matrix = gsl_matrix_calloc(rows, columns);
31};
32
33/** Copy constructor of class GSLMatrix.
34 * Allocates GSL structures and copies components from \a *src.
35 * \param *src source matrix
36 */
37GSLMatrix::GSLMatrix(const GSLMatrix * const src) : rows(src->rows), columns(src->columns)
38{
39 matrix = gsl_matrix_alloc(rows, columns);
40 gsl_matrix_memcpy (matrix, src->matrix);
41};
42
43/** Destructor of class GSLMatrix.
44 * Frees GSL structures
45 */
46GSLMatrix::~GSLMatrix()
47{
48 gsl_matrix_free(matrix);
49 rows = 0;
50 columns = 0;
51};
52
53/** Assignment operator.
54 * \param &rhs right hand side
55 * \return object itself
56 */
57GSLMatrix& GSLMatrix::operator=(const GSLMatrix& rhs)
58{
59 if (this == &rhs) // not necessary here, but identity assignment check is generally a good idea
60 return *this;
61 assert(rows == rhs.rows && columns == rhs.columns && "Number of rows and columns do not match!");
62
63 gsl_matrix_memcpy (matrix, rhs.matrix);
64
65 return *this;
66};
67
68/* ============================ Accessing =============================== */
69/** This function sets the matrix from a double array.
70 * Creates a matrix view of the array and performs a memcopy.
71 * \param *x array of values (no dimension check is performed)
72 */
73void GSLMatrix::SetFromDoubleArray(double * x)
74{
75 gsl_matrix_view m = gsl_matrix_view_array (x, rows, columns);
76 gsl_matrix_memcpy (matrix, &m.matrix);
77};
78
79/** This function returns the i-th element of a matrix.
80 * If \a m or \a n lies outside the allowed range of 0 to GSLMatrix::dimension-1 then the error handler is invoked and 0 is returned.
81 * \param m row index
82 * \param n colum index
83 * \return (m,n)-th element of matrix
84 */
85double GSLMatrix::Get(size_t m, size_t n)
86{
87 return gsl_matrix_get (matrix, m, n);
88};
89
90/** This function sets the value of the \a m -th element of a matrix to \a x.
91 * If \a m or \a n lies outside the allowed range of 0 to GSLMatrix::dimension-1 then the error handler is invoked.
92 * \param m row index
93 * \param m column index
94 * \param x value to set element (m,n) to
95 */
96void GSLMatrix::Set(size_t m, size_t n, double x)
97{
98 gsl_matrix_set (matrix, m, n, x);
99};
100
101/** These functions return a pointer to the \a m-th element of a matrix.
102 * If \a m or \a n lies outside the allowed range of 0 to GSLMatrix::dimension-1 then the error handler is invoked and a null pointer is returned.
103 * \param m index
104 * \return pointer to \a m-th element
105 */
106double *GSLMatrix::Pointer(size_t m, size_t n)
107{
108 return gsl_matrix_ptr (matrix, m, n);
109};
110
111/** These functions return a constant pointer to the \a m-th element of a matrix.
112 * If \a m or \a n lies outside the allowed range of 0 to GSLMatrix::dimension-1 then the error handler is invoked and a null pointer is returned.
113 * \param m index
114 * \return const pointer to \a m-th element
115 */
116const double *GSLMatrix::const_Pointer(size_t m, size_t n)
117{
118 return gsl_matrix_const_ptr (matrix, m, n);
119};
120
121/* ========================== Initializing =============================== */
122/** This function sets all the elements of the matrix to the value \a x.
123 * \param *x
124 */
125void GSLMatrix::SetAll(double x)
126{
127 gsl_matrix_set_all (matrix, x);
128};
129
130/** This function sets all the elements of the matrix to zero.
131 */
132void GSLMatrix::SetZero()
133{
134 gsl_matrix_set_zero (matrix);
135};
136
137/** This function sets the elements of the matrix to the corresponding elements of the identity matrix, \f$m(i,j) = \delta(i,j)\f$, i.e. a unit diagonal with all off-diagonal elements zero.
138 * This applies to both square and rectangular matrices.
139 */
140void GSLMatrix::SetIdentity()
141{
142 gsl_matrix_set_identity (matrix);
143};
144
145/* ====================== Exchanging elements ============================ */
146/** This function exchanges the \a i-th and \a j-th row of the matrix in-place.
147 * \param i i-th row to swap with ...
148 * \param j ... j-th row to swap against
149 */
150bool GSLMatrix::SwapRows(size_t i, size_t j)
151{
152 return (gsl_matrix_swap_rows (matrix, i, j) == GSL_SUCCESS);
153};
154
155/** This function exchanges the \a i-th and \a j-th column of the matrix in-place.
156 * \param i i-th column to swap with ...
157 * \param j ... j-th column to swap against
158 */
159bool GSLMatrix::SwapColumns(size_t i, size_t j)
160{
161 return (gsl_matrix_swap_columns (matrix, i, j) == GSL_SUCCESS);
162};
163
164/** This function exchanges the \a i-th row and \a j-th column of the matrix in-place.
165 * The matrix must be square for this operation to be possible.
166 * \param i i-th row to swap with ...
167 * \param j ... j-th column to swap against
168 */
169bool GSLMatrix::SwapRowColumn(size_t i, size_t j)
170{
171 assert (rows == columns && "The matrix must be square for swapping row against column to be possible.");
172 return (gsl_matrix_swap_rowcol (matrix, i, j) == GSL_SUCCESS);
173};
174
175/** This function transposes the matrix.
176 * Note that the function is extended to the non-square case, where the matrix is re-allocated and copied.
177 */
178bool GSLMatrix::Transpose()
179{
180 if (rows == columns)// if square, use GSL
181 return (gsl_matrix_transpose (matrix) == GSL_SUCCESS);
182 else { // otherwise we have to copy a bit
183 gsl_matrix *dest = gsl_matrix_alloc(columns, rows);
184 for (size_t i=0;i<rows; i++)
185 for (size_t j=0;j<columns;j++) {
186 gsl_matrix_set(dest, j,i, gsl_matrix_get(matrix, i,j) );
187 }
188 gsl_matrix_free(matrix);
189 matrix = dest;
190 flip(rows, columns);
191 return true;
192 }
193};
194
195/* ============================ Properties ============================== */
196/** Checks whether matrix' elements are strictly null.
197 * \return true - is null, false - else
198 */
199bool GSLMatrix::IsNull()
200{
201 return gsl_matrix_isnull (matrix);
202};
203
204/** Checks whether matrix' elements are strictly positive.
205 * \return true - is positive, false - else
206 */
207bool GSLMatrix::IsPositive()
208{
209 return gsl_matrix_ispos (matrix);
210};
211
212/** Checks whether matrix' elements are strictly negative.
213 * \return true - is negative, false - else
214 */
215bool GSLMatrix::IsNegative()
216{
217 return gsl_matrix_isneg (matrix);
218};
219
220/** Checks whether matrix' elements are strictly non-negative.
221 * \return true - is non-negative, false - else
222 */
223bool GSLMatrix::IsNonNegative()
224{
225 return gsl_matrix_isnonneg (matrix);
226};
227
228/** This function performs a Cholesky decomposition to determine whether matrix is positive definite.
229 * We check whether GSL returns GSL_EDOM as error, indicating that decomposition failed due to matrix not being positive-definite.
230 * \return true - matrix is positive-definite, false - else
231 */
232bool GSLMatrix::IsPositiveDefinite()
233{
234 if (rows != columns) // only possible for square matrices.
235 return false;
236 else
237 return (gsl_linalg_cholesky_decomp (matrix) != GSL_EDOM);
238};
239
240
241/** Calculates the determinant of the matrix.
242 * if matrix is square, uses LU decomposition.
243 */
244double GSLMatrix::Determinant()
245{
246 int signum = 0;
247 assert (rows == columns && "Determinant can only be calculated for square matrices.");
248 gsl_permutation *p = gsl_permutation_alloc(rows);
249 gsl_linalg_LU_decomp(matrix, p, &signum);
250 gsl_permutation_free(p);
251 return gsl_linalg_LU_det(matrix, signum);
252};
253
Note: See TracBrowser for help on using the repository browser.