Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2013-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 "core/ActionShortcut.h" 23 : #include "core/ActionRegister.h" 24 : #include "MultiColvarShortcuts.h" 25 : 26 : #include <string> 27 : #include <cmath> 28 : 29 : namespace PLMD { 30 : namespace multicolvar { 31 : 32 : //+PLUMEDOC COLVAR ANGLES 33 : /* 34 : Calculate an angle. 35 : 36 : This command can be used to compute the angle between three atoms. Alternatively 37 : if four atoms appear in the atom 38 : specification it calculates the angle between 39 : two vectors identified by two pairs of atoms. 40 : 41 : If _three_ atoms are given, the angle is defined as: 42 : \f[ 43 : \theta=\arccos\left(\frac{ {\bf r}_{21}\cdot {\bf r}_{23}}{ 44 : |{\bf r}_{21}| |{\bf r}_{23}|}\right) 45 : \f] 46 : Here \f$ {\bf r}_{ij}\f$ is the distance vector among the 47 : i-th and the j-th listed atom. 48 : 49 : If _four_ atoms are given, the angle is defined as: 50 : \f[ 51 : \theta=\arccos\left(\frac{ {\bf r}_{21}\cdot {\bf r}_{34}}{ 52 : |{\bf r}_{21}| |{\bf r}_{34}|}\right) 53 : \f] 54 : 55 : Notice that angles defined in this way are non-periodic variables and 56 : their value is limited by definition between 0 and \f$\pi\f$. 57 : 58 : The vectors \f$ {\bf r}_{ij}\f$ are by default evaluated taking 59 : periodic boundary conditions into account. 60 : This behavior can be changed with the NOPBC flag. 61 : 62 : \par Examples 63 : 64 : This command tells plumed to calculate the angle between the vector connecting atom 1 to atom 2 and 65 : the vector connecting atom 2 to atom 3 and to print it on file COLVAR1. At the same time, 66 : the angle between vector connecting atom 1 to atom 2 and the vector connecting atom 3 to atom 4 is printed 67 : on file COLVAR2. 68 : \plumedfile 69 : 70 : a: ANGLE ATOMS=1,2,3 71 : # equivalently one could state: 72 : # a: ANGLE ATOMS=1,2,2,3 73 : 74 : b: ANGLE ATOMS=1,2,3,4 75 : 76 : PRINT ARG=a FILE=COLVAR1 77 : PRINT ARG=b FILE=COLVAR2 78 : \endplumedfile 79 : 80 : This final example instructs plumed to calculate all the angles in the first coordination 81 : spheres of the atoms. The bins for a normalized histogram of the distribution is then output 82 : 83 : \plumedfile 84 : ANGLES GROUP=1-38 HISTOGRAM={GAUSSIAN LOWER=0.0 UPPER=pi NBINS=20} SWITCH={GAUSSIAN R_0=1.0} LABEL=a1 85 : PRINT ARG=a1.* FILE=colvar 86 : \endplumedfile 87 : 88 : 89 : */ 90 : //+ENDPLUMEDOC 91 : 92 : class Angles : public ActionShortcut { 93 : public: 94 : explicit Angles(const ActionOptions&); 95 : static void registerKeywords( Keywords& keys ); 96 : }; 97 : 98 : PLUMED_REGISTER_ACTION(Angles,"ANGLES") 99 : 100 4 : void Angles::registerKeywords( Keywords& keys ) { 101 4 : ActionShortcut::registerKeywords( keys ); 102 8 : keys.add("atoms-1","GROUP","Calculate angles for each distinct set of three atoms in the group"); 103 8 : keys.add("atoms-2","GROUPA","A group of central atoms about which angles should be calculated"); 104 8 : keys.add("atoms-2","GROUPB","When used in conjunction with GROUPA this keyword instructs plumed " 105 : "to calculate all distinct angles involving one atom from GROUPA " 106 : "and two atoms from GROUPB. The atom from GROUPA is the central atom."); 107 8 : keys.add("atoms-3","GROUPC","This must be used in conjunction with GROUPA and GROUPB. All angles " 108 : "involving one atom from GROUPA, one atom from GROUPB and one atom from " 109 : "GROUPC are calculated. The GROUPA atoms are assumed to be the central " 110 : "atoms"); 111 8 : keys.add("optional","SWITCH","the switching function specifies that only those bonds that have a length that is less than a certain threshold are considered"); 112 8 : keys.setValueDescription("vector","the ANGLE for each set of three atoms that were specified"); 113 8 : MultiColvarShortcuts::shortcutKeywords( keys ); keys.needsAction("ANGLE"); keys.needsAction("COORD_ANGLES"); 114 4 : } 115 : 116 1 : Angles::Angles(const ActionOptions&ao): 117 : Action(ao), 118 1 : ActionShortcut(ao) 119 : { 120 2 : std::string swit; parse("SWITCH",swit); 121 1 : if( swit.length()>0 ) { 122 0 : std::string cat, grp; parse("GROUPA",cat); parse("GROUPB",grp); 123 0 : if( cat.length()==0 || grp.length()==0 ) error("must use GROUPA/GROUPB when using SWITCH"); 124 0 : readInputLine( getShortcutLabel() + ": COORD_ANGLES SWITCH={" + swit + "} CATOMS=" + cat + " GROUP=" + grp + " " + convertInputLineToString() ); 125 : return; 126 : } 127 2 : std::vector<std::string> group; parseVector("GROUP",group); 128 2 : std::vector<std::string> groupa; parseVector("GROUPA",groupa); 129 2 : std::vector<std::string> groupb; parseVector("GROUPB",groupb); 130 2 : std::vector<std::string> groupc; parseVector("GROUPC",groupc); 131 1 : if( group.size()>0 ) { 132 0 : if( groupa.size()>0 || groupb.size()>0 || groupc.size()>0 ) error("should only be GROUP keyword in input not GROUPA/GROUPB/GROUPC"); 133 0 : Tools::interpretRanges( group ); std::string ainput = getShortcutLabel() + ": ANGLE"; unsigned n=1; 134 : // Not sure if this triple sum makes any sense 135 0 : for(unsigned i=2; i<group.size(); ++i ) { 136 0 : for(unsigned j=1; j<i; ++j ) { 137 0 : for(unsigned k=0; k<j; ++k) { 138 0 : std::string str_n; Tools::convert( n, str_n ); 139 0 : ainput += " ATOMS" + str_n + "=" + group[i] + "," + group[j] + "," + group[k]; n++; 140 : } 141 : } 142 : } 143 0 : readInputLine( ainput ); 144 1 : } else if( groupc.size()>0 ) { 145 0 : Tools::interpretRanges( groupa ); Tools::interpretRanges( groupb ); Tools::interpretRanges( groupc ); 146 0 : unsigned n=1; std::string ainput = getShortcutLabel() + ": ANGLE"; 147 0 : for(unsigned i=0; i<groupa.size(); ++i ) { 148 0 : for(unsigned j=0; j<groupb.size(); ++j ) { 149 0 : for(unsigned k=0; k<groupc.size(); ++k) { 150 0 : std::string str_n; Tools::convert( n, str_n ); 151 0 : ainput += " ATOMS" + str_n + "=" + groupb[j] + "," + groupa[i] + "," + groupc[k]; n++; 152 : } 153 : } 154 : } 155 0 : readInputLine( ainput ); 156 1 : } else if( groupa.size()>0 ) { 157 1 : Tools::interpretRanges( groupa ); Tools::interpretRanges( groupb ); 158 1 : unsigned n=1; std::string ainput; ainput = getShortcutLabel() + ": ANGLE"; 159 2 : for(unsigned i=0; i<groupa.size(); ++i ) { 160 9 : for(unsigned j=1; j<groupb.size(); ++j ) { 161 44 : for(unsigned k=0; k<j; ++k) { 162 36 : std::string str_n; Tools::convert( n, str_n ); 163 72 : ainput += " ATOMS" + str_n + "=" + groupb[j] + "," + groupa[i] + "," + groupb[k]; n++; 164 : } 165 : } 166 : } 167 1 : readInputLine( ainput ); 168 : } 169 1 : MultiColvarShortcuts::expandFunctions( getShortcutLabel(), getShortcutLabel(), "", this ); 170 1 : } 171 : 172 : } 173 : } 174 : 175 : 176 :