Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2013-2020 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 "core/ActionShortcut.h" 23 : #include "core/ActionRegister.h" 24 : #include "core/PlumedMain.h" 25 : #include "core/ActionSet.h" 26 : #include "core/Group.h" 27 : #include "AdjacencyMatrixBase.h" 28 : 29 : //+PLUMEDOC MATRIX CONTACT_MATRIX 30 : /* 31 : Adjacency matrix in which two atoms are adjacent if they are within a certain cutoff. 32 : 33 : As discussed in the section of the manual on \ref contactmatrix a useful tool for developing complex collective variables is the notion of the 34 : so called adjacency matrix. An adjacency matrix is an \f$N \times N\f$ matrix in which the \f$i\f$th, \f$j\f$th element tells you whether 35 : or not the \f$i\f$th and \f$j\f$th atoms/molecules from a set of \f$N\f$ atoms/molecules are adjacent or not. These matrices can then be further 36 : analyzed using a number of other algorithms as is detailed in \cite tribello-clustering. 37 : 38 : For this action the elements of the contact matrix are calculated using: 39 : 40 : \f[ 41 : a_{ij} = \sigma( |\mathbf{r}_{ij}| ) 42 : \f] 43 : 44 : where \f$|\mathbf{r}_{ij}|\f$ is the magnitude of the vector connecting atoms \f$i\f$ and \f$j\f$ and where \f$\sigma\f$ is a \ref switchingfunction. 45 : 46 : \par Examples 47 : 48 : The input shown below calculates a \f$6 \times 6\f$ matrix whose elements are equal to one if atom \f$i\f$ and atom \f$j\f$ are within 0.3 nm 49 : of each other and which is zero otherwise. The columns in this matrix are then summed so as to give the coordination number for each atom. 50 : The final quantity output in the colvar file is thus the average coordination number. 51 : 52 : \plumedfile 53 : mat: CONTACT_MATRIX ATOMS=1-6 SWITCH={EXP D_0=0.2 R_0=0.1 D_MAX=0.66} 54 : COLUMNSUMS MATRIX=mat MEAN LABEL=csums 55 : PRINT ARG=csums.* FILE=colvar 56 : \endplumedfile 57 : 58 : */ 59 : //+ENDPLUMEDOC 60 : 61 : namespace PLMD { 62 : namespace adjmat { 63 : 64 : class ContactMatrixShortcut : public ActionShortcut { 65 : public: 66 : static void registerKeywords(Keywords& keys); 67 : explicit ContactMatrixShortcut(const ActionOptions&); 68 : }; 69 : 70 : PLUMED_REGISTER_ACTION(ContactMatrixShortcut,"CONTACT_MATRIX") 71 : 72 334 : void ContactMatrixShortcut::registerKeywords(Keywords& keys) { 73 334 : AdjacencyMatrixBase::registerKeywords( keys ); 74 1002 : keys.remove("GROUP"); keys.add("numbered","GROUP","specifies the list of atoms that should be assumed indistinguishable"); 75 668 : keys.add("compulsory","NN","6","The n parameter of the switching function "); 76 668 : keys.add("compulsory","MM","0","The m parameter of the switching function; 0 implies 2*NN"); 77 668 : keys.add("compulsory","D_0","0.0","The d_0 parameter of the switching function"); 78 668 : keys.add("compulsory","R_0","The r_0 parameter of the switching function"); 79 668 : keys.add("numbered","SWITCH","specify the switching function to use between two sets of indistinguishable atoms"); 80 1002 : keys.addActionNameSuffix("_PROPER"); keys.needsAction("TRANSPOSE"); keys.needsAction("CONCATENATE"); 81 334 : } 82 : 83 208 : ContactMatrixShortcut::ContactMatrixShortcut(const ActionOptions& ao): 84 : Action(ao), 85 208 : ActionShortcut(ao) 86 : { 87 208 : std::vector<std::string> grp_str; std::string atomsstr=""; 88 416 : std::vector<std::string> atomsvec; parseVector("ATOMS",atomsvec); 89 208 : if( atomsvec.size()>0 ) { 90 0 : for(unsigned i=0; i<atomsvec.size(); ++i) { 91 0 : Group* gg = plumed.getActionSet().selectWithLabel<Group*>( atomsvec[i] ); 92 0 : if( gg ) grp_str.push_back( atomsvec[i] ); 93 : } 94 0 : if( grp_str.size()!=atomsvec.size() ) { 95 0 : grp_str.resize(0); 96 0 : atomsstr = " ATOMS=" + atomsvec[0]; for(unsigned i=1; i<atomsvec.size(); ++i) atomsstr += "," + atomsvec[i]; 97 : } 98 : } else { 99 : std::string grp_inpt; 100 2 : for(unsigned i=1;; ++i) { 101 420 : if( !parseNumbered("GROUP",i,grp_inpt) ) break; 102 2 : grp_str.push_back( grp_inpt ); 103 : } 104 : } 105 208 : if( grp_str.size()>9 ) error("cannot handle more than 9 groups"); 106 415 : if( grp_str.size()==0 ) { readInputLine( getShortcutLabel() + ": CONTACT_MATRIX_PROPER " + atomsstr + " " + convertInputLineToString() ); return; } 107 : 108 3 : for(unsigned i=0; i<grp_str.size(); ++i) { 109 4 : std::string sw_str, num; Tools::convert( i+1, num ); parseNumbered("SWITCH", (i+1)*10 + 1 + i, sw_str ); 110 2 : if( sw_str.length()==0 ) error("missing SWITCH" + num + num + " keyword"); 111 4 : readInputLine( getShortcutLabel() + num + num + ": CONTACT_MATRIX_PROPER GROUP=" + grp_str[i] + " SWITCH={" + sw_str + "}" ); 112 3 : for(unsigned j=0; j<i; ++j) { 113 2 : std::string sw_str2, jnum; Tools::convert( j+1, jnum ); parseNumbered("SWITCH", (j+1)*10 + 1 + i, sw_str2); 114 1 : if( sw_str2.length()==0 ) error("missing SWITCH" + jnum + num + " keyword"); 115 2 : readInputLine( getShortcutLabel() + jnum + num + ": CONTACT_MATRIX_PROPER GROUPA=" + grp_str[j] + " GROUPB=" + grp_str[i] + " SWITCH={" + sw_str2 +"}"); 116 2 : readInputLine( getShortcutLabel() + num + jnum + ": TRANSPOSE ARG=" + getShortcutLabel() + jnum + num ); 117 : } 118 : } 119 1 : std::string join_matrices = getShortcutLabel() + ": CONCATENATE"; 120 3 : for(unsigned i=0; i<grp_str.size(); ++i) { 121 2 : std::string inum; Tools::convert(i+1,inum); 122 6 : for(unsigned j=0; j<grp_str.size(); ++j) { 123 4 : std::string jnum; Tools::convert(j+1,jnum); 124 5 : if( i>j ) join_matrices += " MATRIX" + inum + jnum + "=" + getShortcutLabel() + inum + jnum; 125 6 : else join_matrices += " MATRIX" + inum + jnum + "=" + getShortcutLabel() + inum + jnum; 126 : } 127 : } 128 1 : readInputLine( join_matrices ); 129 416 : } 130 : 131 : } 132 : }