Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2015-2023 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 "AdjacencyMatrixBase.h" 23 : #include "core/ActionRegister.h" 24 : #include "tools/SwitchingFunction.h" 25 : #include "tools/Matrix.h" 26 : 27 : //+PLUMEDOC MATRIX CONTACT_MATRIX_PROPER 28 : /* 29 : Adjacency matrix in which two atoms are adjacent if they are within a certain cutoff. 30 : 31 : As discussed in the section of the manual on \ref contactmatrix a useful tool for developing complex collective variables is the notion of the 32 : 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 33 : 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 34 : analyzed using a number of other algorithms as is detailed in \cite tribello-clustering. 35 : 36 : For this action the elements of the contact matrix are calculated using: 37 : 38 : \f[ 39 : a_{ij} = \sigma( |\mathbf{r}_{ij}| ) 40 : \f] 41 : 42 : 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. 43 : 44 : \par Examples 45 : 46 : 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 47 : 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. 48 : The final quantity output in the colvar file is thus the average coordination number. 49 : 50 : \plumedfile 51 : mat: CONTACT_MATRIX ATOMS=1-6 SWITCH={EXP D_0=0.2 R_0=0.1 D_MAX=0.66} 52 : COLUMNSUMS MATRIX=mat MEAN LABEL=csums 53 : PRINT ARG=csums.* FILE=colvar 54 : \endplumedfile 55 : 56 : */ 57 : //+ENDPLUMEDOC 58 : 59 : 60 : namespace PLMD { 61 : namespace adjmat { 62 : 63 : class ContactMatrix : public AdjacencyMatrixBase { 64 : private: 65 : /// switching function 66 : SwitchingFunction switchingFunction; 67 : public: 68 : /// Create manual 69 : static void registerKeywords( Keywords& keys ); 70 : /// Constructor 71 : explicit ContactMatrix(const ActionOptions&); 72 : /// This does nothing 73 : double calculateWeight( const Vector& pos1, const Vector& pos2, const unsigned& natoms, MultiValue& myvals ) const override; 74 : /// Override this so we write the graph properly 75 8 : std::string writeInGraph() const override { return "CONTACT_MATRIX"; } 76 : }; 77 : 78 : PLUMED_REGISTER_ACTION(ContactMatrix,"CONTACT_MATRIX_PROPER") 79 : 80 424 : void ContactMatrix::registerKeywords( Keywords& keys ) { 81 424 : AdjacencyMatrixBase::registerKeywords( keys ); keys.setDisplayName("CONTACT_MATRIX"); 82 848 : keys.add("compulsory","NN","6","The n parameter of the switching function "); 83 848 : keys.add("compulsory","MM","0","The m parameter of the switching function; 0 implies 2*NN"); 84 848 : keys.add("compulsory","D_0","0.0","The d_0 parameter of the switching function"); 85 848 : keys.add("compulsory","R_0","The r_0 parameter of the switching function"); 86 848 : keys.add("optional","SWITCH","This keyword is used if you want to employ an alternative to the continuous swiching function defined above. " 87 : "The following provides information on the \\ref switchingfunction that are available. " 88 : "When this keyword is present you no longer need the NN, MM, D_0 and R_0 keywords."); 89 424 : } 90 : 91 210 : ContactMatrix::ContactMatrix( const ActionOptions& ao ): 92 : Action(ao), 93 210 : AdjacencyMatrixBase(ao) 94 : { 95 420 : std::string errors, input; parse("SWITCH",input); 96 210 : if( input.length()>0 ) { 97 183 : switchingFunction.set( input, errors ); 98 183 : if( errors.length()!=0 ) error("problem reading switching function description " + errors); 99 : } else { 100 27 : double r_0=-1.0, d_0; int nn, mm; 101 54 : parse("NN",nn); parse("MM",mm); 102 54 : parse("R_0",r_0); parse("D_0",d_0); 103 27 : if( r_0<0.0 ) error("you must set a value for R_0"); 104 27 : switchingFunction.set(nn,mm,r_0,d_0); 105 : } 106 : // And set the link cell cutoff 107 210 : log.printf(" switching function cutoff is %s \n",switchingFunction.description().c_str() ); 108 210 : setLinkCellCutoff( true, switchingFunction.get_dmax() ); 109 210 : } 110 : 111 15298242 : double ContactMatrix::calculateWeight( const Vector& pos1, const Vector& pos2, const unsigned& natoms, MultiValue& myvals ) const { 112 15298242 : Vector distance = pos2; double mod2 = distance.modulo2(); 113 15298242 : if( mod2<epsilon ) return 0.0; // Atoms can't be bonded to themselves 114 15296728 : double dfunc, val = switchingFunction.calculateSqr( mod2, dfunc ); 115 15296728 : if( val<epsilon ) return 0.0; 116 8264002 : if( doNotCalculateDerivatives() ) return val; 117 5500819 : addAtomDerivatives( 0, (-dfunc)*distance, myvals ); 118 5500819 : addAtomDerivatives( 1, (+dfunc)*distance, myvals ); 119 5500819 : addBoxDerivatives( (-dfunc)*Tensor(distance,distance), myvals ); 120 5500819 : return val; 121 : } 122 : 123 : } 124 : } 125 :