source: src/Actions/MoleculeAction/RotateToPrincipalAxisSystemAction.cpp@ 8f4df1

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

Merge branch 'AtomicPositionEncapsulation' into stable

Conflicts:

src/Actions/AtomAction/ChangeElementAction.cpp
src/Actions/WorldAction/RemoveSphereOfAtomsAction.cpp
src/Makefile.am
src/UIElements/TextUI/TextDialog.cpp
src/analysis_correlation.hpp
src/atom.cpp
src/atom_atominfo.hpp
src/bond.cpp
src/boundary.cpp
src/molecule_geometry.cpp
src/tesselation.cpp
src/tesselationhelpers.cpp
src/triangleintersectionlist.cpp
src/unittests/Makefile.am

  • fixed #includes due to moves to Helpers and LinearAlgebra
  • moved VectorInterface.* and vector_ops.* to LinearAlgebra
  • no more direct access of atom::node, remapped to set/getPosition()
  • no more direct access to atom::type, remapped to set/getType() (also in atom due to derivation and atominfo::AtomicElement is private not protected).
  • Property mode set to 100644
File size: 6.2 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)->getPosition();
86 x -= *CenterOfGravity;
87 const double mass = (*iter)->getType()->mass;
88 InertiaTensor.at(0,0) += mass*(x[1]*x[1] + x[2]*x[2]);
89 InertiaTensor.at(0,1) += mass*(-x[0]*x[1]);
90 InertiaTensor.at(0,2) += mass*(-x[0]*x[2]);
91 InertiaTensor.at(1,0) += mass*(-x[1]*x[0]);
92 InertiaTensor.at(1,1) += mass*(x[0]*x[0] + x[2]*x[2]);
93 InertiaTensor.at(1,2) += mass*(-x[1]*x[2]);
94 InertiaTensor.at(2,0) += mass*(-x[2]*x[0]);
95 InertiaTensor.at(2,1) += mass*(-x[2]*x[1]);
96 InertiaTensor.at(2,2) += mass*(x[0]*x[0] + x[1]*x[1]);
97 }
98 // print InertiaTensor for debugging
99 DoLog(0) && (Log() << Verbose(0) << "The inertia tensor is:" << InertiaTensor << endl);
100
101 // diagonalize to determine principal axis system
102 Vector Eigenvalues = InertiaTensor.transformToEigenbasis();
103
104 for(int i=0;i<NDIM;i++)
105 DoLog(0) && (Log() << Verbose(0) << "eigenvalue = " << Eigenvalues[i] << ", eigenvector = " << InertiaTensor.column(i) << endl);
106
107 // check whether we rotate or not
108 DoLog(0) && (Log() << Verbose(0) << "Transforming molecule into PAS ... ");
109
110 // obtain first column, eigenvector to biggest eigenvalue
111 Vector BiggestEigenvector(InertiaTensor.column(Eigenvalues.SmallestComponent()));
112 Vector DesiredAxis(Axis);
113
114 // Creation Line that is the rotation axis
115 DesiredAxis.VectorProduct(BiggestEigenvector);
116 Line RotationAxis(Vector(0.,0.,0.), DesiredAxis);
117
118 // determine angle
119 const double alpha = BiggestEigenvector.Angle(Axis);
120
121 DoLog(0) && (Log() << Verbose(0) << alpha << endl);
122
123 for (molecule::iterator iter = mol->begin(); iter != mol->end(); ++iter) {
124 *(*iter) -= *CenterOfGravity;
125 (*iter)->setPosition(RotationAxis.rotateVector((*iter)->getPosition(), alpha));
126 *(*iter) += *CenterOfGravity;
127 }
128 DoLog(0) && (Log() << Verbose(0) << "done." << endl);
129
130 // summing anew for debugging (resulting matrix has to be diagonal!)
131 // reset inertia tensor
132 InertiaTensor.zero();
133
134 // sum up inertia tensor
135 for (molecule::const_iterator iter = mol->begin(); iter != mol->end(); ++iter) {
136 Vector x = (*iter)->getPosition();
137 x -= *CenterOfGravity;
138 const double mass = (*iter)->getType()->mass;
139 InertiaTensor.at(0,0) += mass*(x[1]*x[1] + x[2]*x[2]);
140 InertiaTensor.at(0,1) += mass*(-x[0]*x[1]);
141 InertiaTensor.at(0,2) += mass*(-x[0]*x[2]);
142 InertiaTensor.at(1,0) += mass*(-x[1]*x[0]);
143 InertiaTensor.at(1,1) += mass*(x[0]*x[0] + x[2]*x[2]);
144 InertiaTensor.at(1,2) += mass*(-x[1]*x[2]);
145 InertiaTensor.at(2,0) += mass*(-x[2]*x[0]);
146 InertiaTensor.at(2,1) += mass*(-x[2]*x[1]);
147 InertiaTensor.at(2,2) += mass*(x[0]*x[0] + x[1]*x[1]);
148 // print InertiaTensor for debugging
149 DoLog(0) && (Log() << Verbose(0) << "The inertia tensor is:" << InertiaTensor << endl);
150 }
151
152 // free everything
153 delete(CenterOfGravity);
154 }
155 return Action::success;
156}
157
158Action::state_ptr MoleculeRotateToPrincipalAxisSystemAction::performUndo(Action::state_ptr _state) {
159// MoleculeRotateToPrincipalAxisSystemState *state = assert_cast<MoleculeRotateToPrincipalAxisSystemState*>(_state.get());
160
161// string newName = state->mol->getName();
162// state->mol->setName(state->lastName);
163
164 return Action::failure;
165}
166
167Action::state_ptr MoleculeRotateToPrincipalAxisSystemAction::performRedo(Action::state_ptr _state){
168 // Undo and redo have to do the same for this action
169 return performUndo(_state);
170}
171
172bool MoleculeRotateToPrincipalAxisSystemAction::canUndo() {
173 return false;
174}
175
176bool MoleculeRotateToPrincipalAxisSystemAction::shouldUndo() {
177 return false;
178}
179
180const string MoleculeRotateToPrincipalAxisSystemAction::getName() {
181 return NAME;
182}
Note: See TracBrowser for help on using the repository browser.