[bcf653] | 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 |
|
---|
[97ebf8] | 8 | /*
|
---|
| 9 | * RotateToPrincipalAxisSystemAction.cpp
|
---|
| 10 | *
|
---|
| 11 | * Created on: May 10, 2010
|
---|
| 12 | * Author: heber
|
---|
| 13 | */
|
---|
| 14 |
|
---|
[bf3817] | 15 | // include config.h
|
---|
| 16 | #ifdef HAVE_CONFIG_H
|
---|
| 17 | #include <config.h>
|
---|
| 18 | #endif
|
---|
| 19 |
|
---|
[ad011c] | 20 | #include "CodePatterns/MemDebug.hpp"
|
---|
[112b09] | 21 |
|
---|
[ad011c] | 22 | #include "CodePatterns/Log.hpp"
|
---|
| 23 | #include "CodePatterns/Verbose.hpp"
|
---|
[6e5084] | 24 | #include "LinearAlgebra/Line.hpp"
|
---|
[cca9ef] | 25 | #include "LinearAlgebra/RealSpaceMatrix.hpp"
|
---|
[6e5084] | 26 | #include "LinearAlgebra/Vector.hpp"
|
---|
| 27 | #include "element.hpp"
|
---|
| 28 | #include "molecule.hpp"
|
---|
[1a3c26] | 29 |
|
---|
[97ebf8] | 30 |
|
---|
| 31 | #include <iostream>
|
---|
| 32 | #include <fstream>
|
---|
| 33 | #include <string>
|
---|
| 34 |
|
---|
| 35 | using namespace std;
|
---|
| 36 |
|
---|
[1fd675] | 37 | #include "Actions/MoleculeAction/RotateToPrincipalAxisSystemAction.hpp"
|
---|
[845613] | 38 |
|
---|
[1fd675] | 39 | // and construct the stuff
|
---|
| 40 | #include "RotateToPrincipalAxisSystemAction.def"
|
---|
| 41 | #include "Action_impl_pre.hpp"
|
---|
| 42 | /** =========== define the function ====================== */
|
---|
[845613] | 43 | Action::state_ptr MoleculeRotateToPrincipalAxisSystemAction::performCall() {
|
---|
| 44 | molecule *mol = NULL;
|
---|
[6e5084] | 45 |
|
---|
[1fd675] | 46 | // obtain information
|
---|
| 47 | getParametersfromValueStorage();
|
---|
[97ebf8] | 48 |
|
---|
[845613] | 49 | for (World::MoleculeSelectionIterator iter = World::getInstance().beginMoleculeSelection(); iter != World::getInstance().endMoleculeSelection(); ++iter) {
|
---|
| 50 | mol = iter->second;
|
---|
[97ebf8] | 51 | DoLog(0) && (Log() << Verbose(0) << "Converting to prinicipal axis system." << endl);
|
---|
[cca9ef] | 52 | RealSpaceMatrix InertiaTensor;
|
---|
[6e5084] | 53 | Vector *CenterOfGravity = mol->DetermineCenterOfGravity();
|
---|
| 54 |
|
---|
| 55 | // reset inertia tensor
|
---|
[9eb7580] | 56 | InertiaTensor.setZero();
|
---|
[6e5084] | 57 |
|
---|
| 58 | // sum up inertia tensor
|
---|
| 59 | for (molecule::const_iterator iter = mol->begin(); iter != mol->end(); ++iter) {
|
---|
[8f4df1] | 60 | Vector x = (*iter)->getPosition();
|
---|
[6e5084] | 61 | x -= *CenterOfGravity;
|
---|
[83f176] | 62 | const double mass = (*iter)->getType()->getMass();
|
---|
[8f4df1] | 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]);
|
---|
[6e5084] | 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()));
|
---|
[1fd675] | 87 | Vector DesiredAxis(params.Axis);
|
---|
[6e5084] | 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
|
---|
[1fd675] | 94 | const double alpha = BiggestEigenvector.Angle(params.Axis);
|
---|
[6e5084] | 95 |
|
---|
| 96 | DoLog(0) && (Log() << Verbose(0) << alpha << endl);
|
---|
| 97 |
|
---|
| 98 | for (molecule::iterator iter = mol->begin(); iter != mol->end(); ++iter) {
|
---|
[8f4df1] | 99 | *(*iter) -= *CenterOfGravity;
|
---|
| 100 | (*iter)->setPosition(RotationAxis.rotateVector((*iter)->getPosition(), alpha));
|
---|
| 101 | *(*iter) += *CenterOfGravity;
|
---|
[6e5084] | 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
|
---|
[9eb7580] | 107 | InertiaTensor.setZero();
|
---|
[6e5084] | 108 |
|
---|
| 109 | // sum up inertia tensor
|
---|
| 110 | for (molecule::const_iterator iter = mol->begin(); iter != mol->end(); ++iter) {
|
---|
[8f4df1] | 111 | Vector x = (*iter)->getPosition();
|
---|
[6e5084] | 112 | x -= *CenterOfGravity;
|
---|
[83f176] | 113 | const double mass = (*iter)->getType()->getMass();
|
---|
[8f4df1] | 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]);
|
---|
[6e5084] | 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);
|
---|
[97ebf8] | 129 | }
|
---|
[845613] | 130 | return Action::success;
|
---|
[97ebf8] | 131 | }
|
---|
| 132 |
|
---|
| 133 | Action::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 |
|
---|
| 142 | Action::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 |
|
---|
| 147 | bool MoleculeRotateToPrincipalAxisSystemAction::canUndo() {
|
---|
| 148 | return false;
|
---|
| 149 | }
|
---|
| 150 |
|
---|
| 151 | bool MoleculeRotateToPrincipalAxisSystemAction::shouldUndo() {
|
---|
| 152 | return false;
|
---|
| 153 | }
|
---|
[1fd675] | 154 | /** =========== end of function ====================== */
|
---|