Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2014-2017 The plumed team 3 : (see the PEOPLE file at the root of the distribution for a list of names) 4 : 5 : See http://www.plumed.org for more information. 6 : 7 : This file is part of plumed, version 2. 8 : 9 : plumed is free software: you can redistribute it and/or modify 10 : it under the terms of the GNU Lesser General Public License as published by 11 : the Free Software Foundation, either version 3 of the License, or 12 : (at your option) any later version. 13 : 14 : plumed is distributed in the hope that it will be useful, 15 : but WITHOUT ANY WARRANTY; without even the implied warranty of 16 : MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 17 : GNU Lesser General Public License for more details. 18 : 19 : You should have received a copy of the GNU Lesser General Public License 20 : along with plumed. If not, see <http://www.gnu.org/licenses/>. 21 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */ 22 : #include "MatrixOperationBase.h" 23 : #include "core/ActionSetup.h" 24 : #include "core/ActionRegister.h" 25 : 26 : //+PLUMEDOC MCOLVAR INVERT_MATRIX 27 : /* 28 : Calculate the inverse of the input matrix 29 : 30 : This action allows you to calculate the inverse of a real symmetric matrix. 31 : If the matrix is not symmetric then the [dgetrf](https://www.netlib.org/lapack/explore-html-3.6.1/dd/d9a/group__double_g_ecomputational_ga0019443faea08275ca60a734d0593e60.html) 32 : and [dgetri](https://www.netlib.org/lapack/explore-html-3.6.1/dd/d9a/group__double_g_ecomputational_ga56d9c860ce4ce42ded7f914fdb0683ff.html) 33 : functions from the [LAPACK](https://www.netlib.org/lapack/explore-html/) library are used. 34 : If the matrix is symmetric then we use [dsyevr](https://www.netlib.org/lapack/explore-html/d1/d56/group__heevr_gaa334ac0c11113576db0fc37b7565e8b5.html#gaa334ac0c11113576db0fc37b7565e8b5) 35 : to find the eigenvalues. The inverse matrix of the input matrix $M$ with [eigendecomposition](https://en.wikipedia.org/wiki/Eigendecomposition_of_a_matrix).: 36 : 37 : $$ 38 : M = Q \Lambda Q^{-1} 39 : $$ 40 : 41 : is then found as: 42 : 43 : $$ 44 : M^{-1} = Q \Lambda^{-1} Q^{-1} 45 : $$ 46 : 47 : where $\Lambda^{-1}$ is the inverse of the (diagonal) matrix of eigenvalues $\Lambda$ and $Q$ is the matrix of eigenvectors. 48 : 49 : The following example shows how this action is used in practise: 50 : 51 : ```plumed 52 : c: DISTANCE_MATRIX ATOMS=1-4 53 : ci: INVERT_MATRIX ARG=c 54 : PRINT ARG=ci FILE=colvar 55 : ``` 56 : 57 : */ 58 : //+ENDPLUMEDOC 59 : 60 : namespace PLMD { 61 : namespace matrixtools { 62 : 63 : class InvertMatrix : public MatrixOperationBase { 64 : private: 65 : bool input_is_constant; 66 : Matrix<double> mymatrix; 67 : Matrix<double> inverse; 68 : public: 69 : static void registerKeywords( Keywords& keys ); 70 : /// Constructor 71 : explicit InvertMatrix(const ActionOptions&); 72 : /// 73 6 : unsigned getNumberOfDerivatives() override { 74 6 : return 0; 75 : } 76 : /// Do the calculation 77 : void calculate() override; 78 : /// 79 : void apply() override; 80 0 : double getForceOnMatrixElement( const unsigned& jrow, const unsigned& krow ) const override { 81 0 : plumed_error(); 82 : } 83 : }; 84 : 85 : PLUMED_REGISTER_ACTION(InvertMatrix,"INVERT_MATRIX") 86 : 87 18 : void InvertMatrix::registerKeywords( Keywords& keys ) { 88 18 : MatrixOperationBase::registerKeywords( keys ); 89 36 : keys.setValueDescription("matrix","the inverse of the input matrix"); 90 18 : } 91 : 92 9 : InvertMatrix::InvertMatrix(const ActionOptions& ao): 93 : Action(ao), 94 : MatrixOperationBase(ao), 95 9 : input_is_constant(false) { 96 9 : if( getPntrToArgument(0)->getShape()[0]!=getPntrToArgument(0)->getShape()[1] ) { 97 0 : error("input matrix should be square"); 98 : } 99 : 100 9 : ActionSetup* as = dynamic_cast<ActionSetup*>( getPntrToArgument(0)->getPntrToAction() ); 101 9 : if(as) { 102 0 : input_is_constant=true; 103 : } 104 : 105 9 : std::vector<unsigned> shape(2); 106 9 : shape[0]=shape[1]=getPntrToArgument(0)->getShape()[0]; 107 9 : addValue( shape ); 108 9 : setNotPeriodic(); 109 9 : getPntrToComponent(0)->buildDataStore(); 110 9 : getPntrToComponent(0)->reshapeMatrixStore( shape[1] ); 111 9 : mymatrix.resize( shape[0], shape[1] ); 112 9 : inverse.resize( shape[0], shape[1] ); 113 9 : } 114 : 115 9 : void InvertMatrix::calculate() { 116 : // Retrieve the matrix from input 117 9 : retrieveFullMatrix( mymatrix ); 118 : // Now invert the matrix 119 9 : Invert( mymatrix, inverse ); 120 : // And set the inverse 121 : unsigned k = 0; 122 9 : Value* myval=getPntrToComponent(0); 123 30 : for(unsigned i=0; i<mymatrix.nrows(); ++i) { 124 78 : for(unsigned j=0; j<mymatrix.ncols(); ++j) { 125 57 : myval->set( k, inverse(i,j) ); 126 57 : k++; 127 : } 128 : } 129 : 130 9 : if( !doNotCalculateDerivatives() && !input_is_constant ) { 131 0 : error("derivatives of inverse matrix have not been implemented"); 132 : } 133 9 : } 134 : 135 0 : void InvertMatrix::apply() { 136 0 : if( doNotCalculateDerivatives() || input_is_constant ) { 137 0 : return; 138 : } 139 0 : error("derivatives of inverse matrix have not been implemented"); 140 : } 141 : 142 : } 143 : }