source: src/Actions/MoleculeAction/RotateToPrincipalAxisSystemAction.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@…>, 14 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: 6.3 KB
Line 
1/*
2 * RotateToPrincipalAxisSystemAction.cpp
3 *
4 * Created on: May 10, 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
15#include "Actions/MoleculeAction/RotateToPrincipalAxisSystemAction.hpp"
16#include "Actions/ActionRegistry.hpp"
17#include "Helpers/Log.hpp"
18#include "Helpers/Verbose.hpp"
19#include "LinearAlgebra/Line.hpp"
20#include "LinearAlgebra/Matrix.hpp"
21#include "LinearAlgebra/Vector.hpp"
22#include "element.hpp"
23#include "molecule.hpp"
24
25
26#include <iostream>
27#include <fstream>
28#include <string>
29
30using namespace std;
31
32#include "UIElements/UIFactory.hpp"
33#include "UIElements/Dialog.hpp"
34#include "Actions/ValueStorage.hpp"
35
36/****** MoleculeRotateToPrincipalAxisSystemAction *****/
37
38// memento to remember the state when undoing
39
40//class MoleculeRotateToPrincipalAxisSystemState : public ActionState {
41//public:
42// MoleculeRotateToPrincipalAxisSystemState(molecule* _mol,std::string _lastName) :
43// mol(_mol),
44// lastName(_lastName)
45// {}
46// molecule* mol;
47// std::string lastName;
48//};
49
50const char MoleculeRotateToPrincipalAxisSystemAction::NAME[] = "rotate-to-pas";
51
52MoleculeRotateToPrincipalAxisSystemAction::MoleculeRotateToPrincipalAxisSystemAction() :
53 Action(NAME)
54{}
55
56MoleculeRotateToPrincipalAxisSystemAction::~MoleculeRotateToPrincipalAxisSystemAction()
57{}
58
59void MoleculeRotateToPrincipalAxisSystem(Vector &Axis) {
60 ValueStorage::getInstance().setCurrentValue(MoleculeRotateToPrincipalAxisSystemAction::NAME, Axis);
61 ActionRegistry::getInstance().getActionByName(MoleculeRotateToPrincipalAxisSystemAction::NAME)->call(Action::NonInteractive);
62};
63
64Dialog* MoleculeRotateToPrincipalAxisSystemAction::fillDialog(Dialog *dialog) {
65 ASSERT(dialog,"No Dialog given when filling action dialog");
66
67 dialog->queryVector(NAME, false, ValueStorage::getInstance().getDescription(NAME));
68
69 return dialog;
70}
71
72Action::state_ptr MoleculeRotateToPrincipalAxisSystemAction::performCall() {
73 molecule *mol = NULL;
74 Vector Axis;
75
76 // obtain axis to rotate to
77 ValueStorage::getInstance().queryCurrentValue(NAME, Axis);
78
79 for (World::MoleculeSelectionIterator iter = World::getInstance().beginMoleculeSelection(); iter != World::getInstance().endMoleculeSelection(); ++iter) {
80 mol = iter->second;
81 DoLog(0) && (Log() << Verbose(0) << "Converting to prinicipal axis system." << endl);
82 Matrix InertiaTensor;
83 Vector *CenterOfGravity = mol->DetermineCenterOfGravity();
84
85 // reset inertia tensor
86 InertiaTensor.zero();
87
88 // sum up inertia tensor
89 for (molecule::const_iterator iter = mol->begin(); iter != mol->end(); ++iter) {
90 Vector x = (*iter)->getPosition();
91 x -= *CenterOfGravity;
92 const double mass = (*iter)->getType()->mass;
93 InertiaTensor.at(0,0) += mass*(x[1]*x[1] + x[2]*x[2]);
94 InertiaTensor.at(0,1) += mass*(-x[0]*x[1]);
95 InertiaTensor.at(0,2) += mass*(-x[0]*x[2]);
96 InertiaTensor.at(1,0) += mass*(-x[1]*x[0]);
97 InertiaTensor.at(1,1) += mass*(x[0]*x[0] + x[2]*x[2]);
98 InertiaTensor.at(1,2) += mass*(-x[1]*x[2]);
99 InertiaTensor.at(2,0) += mass*(-x[2]*x[0]);
100 InertiaTensor.at(2,1) += mass*(-x[2]*x[1]);
101 InertiaTensor.at(2,2) += mass*(x[0]*x[0] + x[1]*x[1]);
102 }
103 // print InertiaTensor for debugging
104 DoLog(0) && (Log() << Verbose(0) << "The inertia tensor is:" << InertiaTensor << endl);
105
106 // diagonalize to determine principal axis system
107 Vector Eigenvalues = InertiaTensor.transformToEigenbasis();
108
109 for(int i=0;i<NDIM;i++)
110 DoLog(0) && (Log() << Verbose(0) << "eigenvalue = " << Eigenvalues[i] << ", eigenvector = " << InertiaTensor.column(i) << endl);
111
112 // check whether we rotate or not
113 DoLog(0) && (Log() << Verbose(0) << "Transforming molecule into PAS ... ");
114
115 // obtain first column, eigenvector to biggest eigenvalue
116 Vector BiggestEigenvector(InertiaTensor.column(Eigenvalues.SmallestComponent()));
117 Vector DesiredAxis(Axis);
118
119 // Creation Line that is the rotation axis
120 DesiredAxis.VectorProduct(BiggestEigenvector);
121 Line RotationAxis(Vector(0.,0.,0.), DesiredAxis);
122
123 // determine angle
124 const double alpha = BiggestEigenvector.Angle(Axis);
125
126 DoLog(0) && (Log() << Verbose(0) << alpha << endl);
127
128 for (molecule::iterator iter = mol->begin(); iter != mol->end(); ++iter) {
129 *(*iter) -= *CenterOfGravity;
130 (*iter)->setPosition(RotationAxis.rotateVector((*iter)->getPosition(), alpha));
131 *(*iter) += *CenterOfGravity;
132 }
133 DoLog(0) && (Log() << Verbose(0) << "done." << endl);
134
135 // summing anew for debugging (resulting matrix has to be diagonal!)
136 // reset inertia tensor
137 InertiaTensor.zero();
138
139 // sum up inertia tensor
140 for (molecule::const_iterator iter = mol->begin(); iter != mol->end(); ++iter) {
141 Vector x = (*iter)->getPosition();
142 x -= *CenterOfGravity;
143 const double mass = (*iter)->getType()->mass;
144 InertiaTensor.at(0,0) += mass*(x[1]*x[1] + x[2]*x[2]);
145 InertiaTensor.at(0,1) += mass*(-x[0]*x[1]);
146 InertiaTensor.at(0,2) += mass*(-x[0]*x[2]);
147 InertiaTensor.at(1,0) += mass*(-x[1]*x[0]);
148 InertiaTensor.at(1,1) += mass*(x[0]*x[0] + x[2]*x[2]);
149 InertiaTensor.at(1,2) += mass*(-x[1]*x[2]);
150 InertiaTensor.at(2,0) += mass*(-x[2]*x[0]);
151 InertiaTensor.at(2,1) += mass*(-x[2]*x[1]);
152 InertiaTensor.at(2,2) += mass*(x[0]*x[0] + x[1]*x[1]);
153 // print InertiaTensor for debugging
154 DoLog(0) && (Log() << Verbose(0) << "The inertia tensor is:" << InertiaTensor << endl);
155 }
156
157 // free everything
158 delete(CenterOfGravity);
159 }
160 return Action::success;
161}
162
163Action::state_ptr MoleculeRotateToPrincipalAxisSystemAction::performUndo(Action::state_ptr _state) {
164// MoleculeRotateToPrincipalAxisSystemState *state = assert_cast<MoleculeRotateToPrincipalAxisSystemState*>(_state.get());
165
166// string newName = state->mol->getName();
167// state->mol->setName(state->lastName);
168
169 return Action::failure;
170}
171
172Action::state_ptr MoleculeRotateToPrincipalAxisSystemAction::performRedo(Action::state_ptr _state){
173 // Undo and redo have to do the same for this action
174 return performUndo(_state);
175}
176
177bool MoleculeRotateToPrincipalAxisSystemAction::canUndo() {
178 return false;
179}
180
181bool MoleculeRotateToPrincipalAxisSystemAction::shouldUndo() {
182 return false;
183}
184
185const string MoleculeRotateToPrincipalAxisSystemAction::getName() {
186 return NAME;
187}
Note: See TracBrowser for help on using the repository browser.