source: src/Actions/MoleculeAction/RotateToPrincipalAxisSystemAction.cpp@ 6e5084

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

Fixed PrincipalAxisSystemAction.

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