source: src/LinearAlgebra/RealSpaceMatrix.cpp@ 66fd49

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

Rewrite of FillVoidWithMoleculeAction.

FillVoidWithMoleculeAction:

  • new parameter MinDistance and default value of 0.
  • BUGFIX: filler is already created when parsing file, removed useless creation of it initially (also caused lots of confusion due to an "extra" molecule).
  • Undo implemented, regression test inserted.
  • Redo is somewhat hard to implement, as one would use performCall() if it only it would not retrieve its values from ValueStorage ...

FillVoidWithMolecule():

  • filler is now the zeroth not the last molecule, marked by firstInsertion and firstInserter. Filler is removed if no molecules are filled.
  • outsourced stuff into smaller functions
  • removed FillIt to through every atom despite only CurrentPosition, indepedent of atom position, is checked.

TESTFIXES:

  • Analysis/3: test.xyz changed because boundary is now 1.5 instead of 2.1 as 2.1 is not enough of molecules get filled in (and the filler already is).
  • Analysis/3: tensid.data was actually lacking water at (0,0,0) which is after the rewrite present.
  • Property mode set to 100644
File size: 9.0 KB
Line 
1/*
2 * Project: MoleCuilder
3 * Description: creates and alters molecular systems
4 * Copyright (C) 2010 University of Bonn. All rights reserved.
5 * Please see the LICENSE file or "Copyright notice" in builder.cpp for details.
6 */
7
8/*
9 * RealSpaceMatrix.cpp
10 *
11 * Created on: Jun 25, 2010
12 * Author: crueger
13 */
14
15// include config.h
16#ifdef HAVE_CONFIG_H
17#include <config.h>
18#endif
19
20#include "CodePatterns/MemDebug.hpp"
21
22#include "Exceptions/NotInvertibleException.hpp"
23#include "CodePatterns/Assert.hpp"
24#include "Helpers/defs.hpp"
25#include "Helpers/fast_functions.hpp"
26#include "LinearAlgebra/MatrixContent.hpp"
27#include "LinearAlgebra/RealSpaceMatrix.hpp"
28#include "LinearAlgebra/Vector.hpp"
29#include "LinearAlgebra/VectorContent.hpp"
30
31#include <gsl/gsl_blas.h>
32#include <gsl/gsl_eigen.h>
33#include <gsl/gsl_matrix.h>
34#include <gsl/gsl_multimin.h>
35#include <gsl/gsl_vector.h>
36#include <cmath>
37#include <iostream>
38
39using namespace std;
40
41RealSpaceMatrix::RealSpaceMatrix()
42{
43 content = new MatrixContent(NDIM, NDIM);
44 createViews();
45}
46
47RealSpaceMatrix::RealSpaceMatrix(const double* src)
48{
49 content = new MatrixContent(NDIM, NDIM, src);
50 createViews();
51}
52
53RealSpaceMatrix::RealSpaceMatrix(const RealSpaceMatrix &src)
54{
55 content = new MatrixContent(src.content);
56 createViews();
57}
58
59RealSpaceMatrix::RealSpaceMatrix(const MatrixContent &src)
60{
61 content = new MatrixContent(src);
62 createViews();
63}
64
65RealSpaceMatrix::RealSpaceMatrix(MatrixContent* _content)
66{
67 content = new MatrixContent(_content);
68 createViews();
69}
70
71RealSpaceMatrix::~RealSpaceMatrix()
72{
73 // delete all views
74 for(int i=NDIM;i--;){
75 delete rows_ptr[i];
76 }
77 for(int i=NDIM;i--;){
78 delete columns_ptr[i];
79 }
80 delete diagonal_ptr;
81 delete content;
82}
83
84void RealSpaceMatrix::createViews(){
85 // create row views
86 for(int i=NDIM;i--;){
87 VectorContent *rowContent = content->getRowVector(i);
88 rows_ptr[i] = new Vector(rowContent);
89 ASSERT(rowContent == NULL, "RealSpaceMatrix::createViews() - rowContent was not taken over as supposed to happen!");
90 }
91 // create column views
92 for(int i=NDIM;i--;){
93 VectorContent *columnContent = content->getColumnVector(i);
94 columns_ptr[i] = new Vector(columnContent);
95 ASSERT(columnContent == NULL, "RealSpaceMatrix::createViews() - columnContent was not taken over as supposed to happen!");
96 }
97 // create diagonal view
98 VectorContent *diagonalContent = content->getDiagonalVector();
99 diagonal_ptr = new Vector(diagonalContent);
100 ASSERT(diagonalContent == NULL, "RealSpaceMatrix::createViews() - diagonalContent was not taken over as supposed to happen!");
101}
102
103void RealSpaceMatrix::setIdentity(){
104 content->setIdentity();
105}
106
107void RealSpaceMatrix::setZero(){
108 content->setZero();
109}
110
111void RealSpaceMatrix::setRotation(const double x, const double y, const double z)
112{
113 set(0,0, cos(y)*cos(z));
114 set(0,1, cos(z)*sin(x)*sin(y) - cos(x)*sin(z));
115 set(0,2, cos(x)*cos(z)*sin(y) + sin(x) * sin(z));
116 set(1,0, cos(y)*sin(z));
117 set(1,1, cos(x)*cos(z) + sin(x)*sin(y)*sin(z));
118 set(1,2, -cos(z)*sin(x) + cos(x)*sin(y)*sin(z));
119 set(2,0, -sin(y));
120 set(2,1, cos(y)*sin(x));
121 set(2,2, cos(x)*cos(y));
122}
123
124void RealSpaceMatrix::setRandomRotation()
125{
126 double phi[NDIM];
127
128 for (int i=0;i<NDIM;i++) {
129 phi[i] = rand()/(RAND_MAX/(2.*M_PI));
130 }
131
132 set(0,0, cos(phi[0]) *cos(phi[2]) + (sin(phi[0])*sin(phi[1])*sin(phi[2])));
133 set(0,1, sin(phi[0]) *cos(phi[2]) - (cos(phi[0])*sin(phi[1])*sin(phi[2])));
134 set(0,2, cos(phi[1])*sin(phi[2]) );
135 set(1,0, -sin(phi[0])*cos(phi[1]) );
136 set(1,1, cos(phi[0])*cos(phi[1]) );
137 set(1,2, sin(phi[1]) );
138 set(2,0, -cos(phi[0]) *sin(phi[2]) + (sin(phi[0])*sin(phi[1])*cos(phi[2])));
139 set(2,1, -sin(phi[0]) *sin(phi[2]) - (cos(phi[0])*sin(phi[1])*cos(phi[2])));
140 set(2,2, cos(phi[1])*cos(phi[2]) );
141}
142
143
144RealSpaceMatrix &RealSpaceMatrix::operator=(const RealSpaceMatrix &src)
145{
146 // MatrixContent checks for self-assignment
147 *content = *(src.content);
148 return *this;
149}
150
151const RealSpaceMatrix &RealSpaceMatrix::operator+=(const RealSpaceMatrix &rhs)
152{
153 *content += *(rhs.content);
154 return *this;
155}
156
157const RealSpaceMatrix &RealSpaceMatrix::operator-=(const RealSpaceMatrix &rhs)
158{
159 *content -= *(rhs.content);
160 return *this;
161}
162
163const RealSpaceMatrix &RealSpaceMatrix::operator*=(const RealSpaceMatrix &rhs)
164{
165 (*content) *= (*rhs.content);
166 return *this;
167}
168
169const RealSpaceMatrix RealSpaceMatrix::operator+(const RealSpaceMatrix &rhs) const
170{
171 RealSpaceMatrix tmp = *this;
172 tmp+=rhs;
173 return tmp;
174}
175
176const RealSpaceMatrix RealSpaceMatrix::operator-(const RealSpaceMatrix &rhs) const
177{
178 RealSpaceMatrix tmp = *this;
179 tmp-=rhs;
180 return tmp;
181}
182
183const RealSpaceMatrix RealSpaceMatrix::operator*(const RealSpaceMatrix &rhs) const
184{
185 RealSpaceMatrix tmp(content);
186 tmp *= rhs;
187 return tmp;
188}
189
190double &RealSpaceMatrix::at(size_t i, size_t j)
191{
192 return content->at(i,j);
193}
194
195const double RealSpaceMatrix::at(size_t i, size_t j) const
196{
197 return content->at(i,j);
198}
199
200Vector &RealSpaceMatrix::row(size_t i)
201{
202 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
203 return *rows_ptr[i];
204}
205
206const Vector &RealSpaceMatrix::row(size_t i) const
207{
208 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
209 return *rows_ptr[i];
210}
211
212Vector &RealSpaceMatrix::column(size_t i)
213{
214 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
215 return *columns_ptr[i];
216}
217
218const Vector &RealSpaceMatrix::column(size_t i) const
219{
220 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
221 return *columns_ptr[i];
222}
223
224Vector &RealSpaceMatrix::diagonal()
225{
226 return *diagonal_ptr;
227}
228
229const Vector &RealSpaceMatrix::diagonal() const
230{
231 return *diagonal_ptr;
232}
233
234void RealSpaceMatrix::set(size_t i, size_t j, const double value)
235{
236 content->set(i,j, value);
237}
238
239double RealSpaceMatrix::determinant() const{
240 return at(0,0)*at(1,1)*at(2,2)
241 + at(0,1)*at(1,2)*at(2,0)
242 + at(0,2)*at(1,0)*at(2,1)
243 - at(2,0)*at(1,1)*at(0,2)
244 - at(2,1)*at(1,2)*at(0,0)
245 - at(2,2)*at(1,0)*at(0,1);
246}
247
248
249RealSpaceMatrix RealSpaceMatrix::invert() const{
250 double det = determinant();
251 if(fabs(det)<MYEPSILON){
252 throw NotInvertibleException(__FILE__,__LINE__);
253 }
254
255 double detReci = 1./det;
256 RealSpaceMatrix res;
257 res.set(0,0, detReci*RDET2(at(1,1),at(2,1),at(1,2),at(2,2))); // A_11
258 res.set(1,0, -detReci*RDET2(at(1,0),at(2,0),at(1,2),at(2,2))); // A_21
259 res.set(2,0, detReci*RDET2(at(1,0),at(2,0),at(1,1),at(2,1))); // A_31
260 res.set(0,1, -detReci*RDET2(at(0,1),at(2,1),at(0,2),at(2,2))); // A_12
261 res.set(1,1, detReci*RDET2(at(0,0),at(2,0),at(0,2),at(2,2))); // A_22
262 res.set(2,1, -detReci*RDET2(at(0,0),at(2,0),at(0,1),at(2,1))); // A_32
263 res.set(0,2, detReci*RDET2(at(0,1),at(1,1),at(0,2),at(1,2))); // A_13
264 res.set(1,2, -detReci*RDET2(at(0,0),at(1,0),at(0,2),at(1,2))); // A_23
265 res.set(2,2, detReci*RDET2(at(0,0),at(1,0),at(0,1),at(1,1))); // A_33
266 return res;
267}
268
269RealSpaceMatrix RealSpaceMatrix::transpose() const
270{
271 RealSpaceMatrix res = RealSpaceMatrix(const_cast<const MatrixContent *>(content)->transpose());
272 return res;
273}
274
275RealSpaceMatrix &RealSpaceMatrix::transpose()
276{
277 content->transpose();
278 return *this;
279}
280
281Vector RealSpaceMatrix::transformToEigenbasis()
282{
283 gsl_vector *eval = content->transformToEigenbasis();
284 Vector evalues(gsl_vector_get(eval,0), gsl_vector_get(eval,1), gsl_vector_get(eval,2));
285 gsl_vector_free(eval);
286 return evalues;
287}
288
289const RealSpaceMatrix &RealSpaceMatrix::operator*=(const double factor)
290 {
291 *content *= factor;
292 return *this;
293}
294
295const RealSpaceMatrix operator*(const double factor,const RealSpaceMatrix& mat)
296{
297 RealSpaceMatrix tmp = mat;
298 return tmp *= factor;
299}
300
301const RealSpaceMatrix operator*(const RealSpaceMatrix &mat,const double factor)
302{
303 return factor*mat;
304}
305
306bool RealSpaceMatrix::operator==(const RealSpaceMatrix &rhs) const
307{
308 return (*content == *(rhs.content));
309}
310
311/** Blows the 6-dimensional \a cell_size array up to a full NDIM by NDIM matrix.
312 * \param *symm 6-dim array of unique symmetric matrix components
313 * \return allocated NDIM*NDIM array with the symmetric matrix
314 */
315RealSpaceMatrix ReturnFullMatrixforSymmetric(const double * const symm)
316{
317 RealSpaceMatrix matrix;
318 matrix.set(0,0, symm[0]);
319 matrix.set(1,0, symm[1]);
320 matrix.set(2,0, symm[3]);
321 matrix.set(0,1, symm[1]);
322 matrix.set(1,1, symm[2]);
323 matrix.set(2,1, symm[4]);
324 matrix.set(0,2, symm[3]);
325 matrix.set(1,2, symm[4]);
326 matrix.set(2,2, symm[5]);
327 return matrix;
328};
329
330ostream &operator<<(ostream &ost,const RealSpaceMatrix &mat)
331{
332 for(int i = 0;i<NDIM;++i){
333 ost << "\n";
334 for(int j = 0; j<NDIM;++j){
335 ost << mat.at(i,j);
336 if(j!=NDIM-1)
337 ost << "; ";
338 }
339 }
340 return ost;
341}
342
343Vector operator*(const RealSpaceMatrix &mat,const Vector &vec)
344{
345 return (*mat.content) * vec;
346}
347
348Vector &operator*=(Vector& lhs,const RealSpaceMatrix &mat)
349{
350 lhs = mat*lhs;
351 return lhs;
352}
353
Note: See TracBrowser for help on using the repository browser.