source: src/Actions/MoleculeAction/RotateToPrincipalAxisSystemAction.cpp@ 0b2ce9

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

Implemented macros for automatically generating repetitive stuff around Actions.

The idea is that only the following items have to be provided for by the user

  • parameters (i.e. for each the following tupel: type, token and reference)
  • perform...(), ...Undo(), ...

Therefore, we have three new files:

  • Action_impl_header.hpp: Is for the declarations in the header
  • Action_impl_pre.hpp: Is before definition of functions
  • Action_impl_post.hpp: is after definition of functions (cleanup)

Changes:

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