source: src/Actions/MoleculeAction/RotateToPrincipalAxisSystemAction.cpp@ 73916f

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

Renamed Matrix to RealSpaceMatrix.

  • class Matrix only represents 3x3 matrices, whereas we are now going to extend the class MatrixContent to contain arbritrary dimensions.
  • renamed class and file
  • Property mode set to 100644
File size: 5.3 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 * RotateToPrincipalAxisSystemAction.cpp
10 *
11 * Created on: May 10, 2010
12 * Author: heber
13 */
14
15// include config.h
16#ifdef HAVE_CONFIG_H
17#include <config.h>
18#endif
19
20#include "Helpers/MemDebug.hpp"
21
22#include "Helpers/Log.hpp"
23#include "Helpers/Verbose.hpp"
24#include "LinearAlgebra/Line.hpp"
25#include "LinearAlgebra/RealSpaceMatrix.hpp"
26#include "LinearAlgebra/Vector.hpp"
27#include "element.hpp"
28#include "molecule.hpp"
29
30
31#include <iostream>
32#include <fstream>
33#include <string>
34
35using namespace std;
36
37#include "Actions/MoleculeAction/RotateToPrincipalAxisSystemAction.hpp"
38
39// and construct the stuff
40#include "RotateToPrincipalAxisSystemAction.def"
41#include "Action_impl_pre.hpp"
42/** =========== define the function ====================== */
43Action::state_ptr MoleculeRotateToPrincipalAxisSystemAction::performCall() {
44 molecule *mol = NULL;
45
46 // obtain information
47 getParametersfromValueStorage();
48
49 for (World::MoleculeSelectionIterator iter = World::getInstance().beginMoleculeSelection(); iter != World::getInstance().endMoleculeSelection(); ++iter) {
50 mol = iter->second;
51 DoLog(0) && (Log() << Verbose(0) << "Converting to prinicipal axis system." << endl);
52 RealSpaceMatrix InertiaTensor;
53 Vector *CenterOfGravity = mol->DetermineCenterOfGravity();
54
55 // reset inertia tensor
56 InertiaTensor.setZero();
57
58 // sum up inertia tensor
59 for (molecule::const_iterator iter = mol->begin(); iter != mol->end(); ++iter) {
60 Vector x = (*iter)->getPosition();
61 x -= *CenterOfGravity;
62 const double mass = (*iter)->getType()->getMass();
63 InertiaTensor.at(0,0) += mass*(x[1]*x[1] + x[2]*x[2]);
64 InertiaTensor.at(0,1) += mass*(-x[0]*x[1]);
65 InertiaTensor.at(0,2) += mass*(-x[0]*x[2]);
66 InertiaTensor.at(1,0) += mass*(-x[1]*x[0]);
67 InertiaTensor.at(1,1) += mass*(x[0]*x[0] + x[2]*x[2]);
68 InertiaTensor.at(1,2) += mass*(-x[1]*x[2]);
69 InertiaTensor.at(2,0) += mass*(-x[2]*x[0]);
70 InertiaTensor.at(2,1) += mass*(-x[2]*x[1]);
71 InertiaTensor.at(2,2) += mass*(x[0]*x[0] + x[1]*x[1]);
72 }
73 // print InertiaTensor for debugging
74 DoLog(0) && (Log() << Verbose(0) << "The inertia tensor is:" << InertiaTensor << endl);
75
76 // diagonalize to determine principal axis system
77 Vector Eigenvalues = InertiaTensor.transformToEigenbasis();
78
79 for(int i=0;i<NDIM;i++)
80 DoLog(0) && (Log() << Verbose(0) << "eigenvalue = " << Eigenvalues[i] << ", eigenvector = " << InertiaTensor.column(i) << endl);
81
82 // check whether we rotate or not
83 DoLog(0) && (Log() << Verbose(0) << "Transforming molecule into PAS ... ");
84
85 // obtain first column, eigenvector to biggest eigenvalue
86 Vector BiggestEigenvector(InertiaTensor.column(Eigenvalues.SmallestComponent()));
87 Vector DesiredAxis(params.Axis);
88
89 // Creation Line that is the rotation axis
90 DesiredAxis.VectorProduct(BiggestEigenvector);
91 Line RotationAxis(Vector(0.,0.,0.), DesiredAxis);
92
93 // determine angle
94 const double alpha = BiggestEigenvector.Angle(params.Axis);
95
96 DoLog(0) && (Log() << Verbose(0) << alpha << endl);
97
98 for (molecule::iterator iter = mol->begin(); iter != mol->end(); ++iter) {
99 *(*iter) -= *CenterOfGravity;
100 (*iter)->setPosition(RotationAxis.rotateVector((*iter)->getPosition(), alpha));
101 *(*iter) += *CenterOfGravity;
102 }
103 DoLog(0) && (Log() << Verbose(0) << "done." << endl);
104
105 // summing anew for debugging (resulting matrix has to be diagonal!)
106 // reset inertia tensor
107 InertiaTensor.setZero();
108
109 // sum up inertia tensor
110 for (molecule::const_iterator iter = mol->begin(); iter != mol->end(); ++iter) {
111 Vector x = (*iter)->getPosition();
112 x -= *CenterOfGravity;
113 const double mass = (*iter)->getType()->getMass();
114 InertiaTensor.at(0,0) += mass*(x[1]*x[1] + x[2]*x[2]);
115 InertiaTensor.at(0,1) += mass*(-x[0]*x[1]);
116 InertiaTensor.at(0,2) += mass*(-x[0]*x[2]);
117 InertiaTensor.at(1,0) += mass*(-x[1]*x[0]);
118 InertiaTensor.at(1,1) += mass*(x[0]*x[0] + x[2]*x[2]);
119 InertiaTensor.at(1,2) += mass*(-x[1]*x[2]);
120 InertiaTensor.at(2,0) += mass*(-x[2]*x[0]);
121 InertiaTensor.at(2,1) += mass*(-x[2]*x[1]);
122 InertiaTensor.at(2,2) += mass*(x[0]*x[0] + x[1]*x[1]);
123 // print InertiaTensor for debugging
124 DoLog(0) && (Log() << Verbose(0) << "The inertia tensor is:" << InertiaTensor << endl);
125 }
126
127 // free everything
128 delete(CenterOfGravity);
129 }
130 return Action::success;
131}
132
133Action::state_ptr MoleculeRotateToPrincipalAxisSystemAction::performUndo(Action::state_ptr _state) {
134// MoleculeRotateToPrincipalAxisSystemState *state = assert_cast<MoleculeRotateToPrincipalAxisSystemState*>(_state.get());
135
136// string newName = state->mol->getName();
137// state->mol->setName(state->lastName);
138
139 return Action::failure;
140}
141
142Action::state_ptr MoleculeRotateToPrincipalAxisSystemAction::performRedo(Action::state_ptr _state){
143 // Undo and redo have to do the same for this action
144 return performUndo(_state);
145}
146
147bool MoleculeRotateToPrincipalAxisSystemAction::canUndo() {
148 return false;
149}
150
151bool MoleculeRotateToPrincipalAxisSystemAction::shouldUndo() {
152 return false;
153}
154/** =========== end of function ====================== */
Note: See TracBrowser for help on using the repository browser.