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 functions of the distribution of angles. 35 : 36 : __This shortcut action allows you to calculate function of a distribution of angles and reproduces the syntax in older PLUMED versions. 37 : If you look at the example inputs below you can 38 : see how the new syntax operates. We would strongly encourage you to use the newer syntax as it offers greater flexibility.__ 39 : 40 : You can use this command to calculate functions, $g$, such as: 41 : 42 : $$ 43 : f(x) = \sum_{ijk} g( \theta_{ijk} ) 44 : $$ 45 : 46 : where $\theta_{ijk}$ is the angle between the vector connecting atom $i$ and and $j$ and the vector connecting atom $j$ and atom $k$. 47 : Alternatively you can use this command to calculate functions such as: 48 : 49 : $$ 50 : f(x) = \sum_{ijk} s(r_{ij})s(r_{jk}) g(\theta_{ijk}) 51 : $$ 52 : 53 : where $s(r)$ is a switching function. This second form means that you can 54 : use this to calculate functions of the angles in the first coordination sphere of 55 : an atom / molecule and hence use the CVs described in the paper in the bibliography below. 56 : 57 : The following example tells plumed to calculate all angles involving 58 : at least one atom from GROUPA and two atoms from GROUPB in which the distances 59 : are less than 1.0. The number of angles between $\frac{\pi}{4}$ and 60 : $\frac{3\pi}{4}$ is then output 61 : 62 : ```plumed 63 : a1: ANGLES GROUPA=1-10 GROUPB=11-100 BETWEEN={GAUSSIAN LOWER=0.25pi UPPER=0.75pi} SWITCH={GAUSSIAN R_0=1.0} 64 : PRINT ARG=a1.between FILE=colvar 65 : ``` 66 : 67 : */ 68 : //+ENDPLUMEDOC 69 : 70 : class Angles : public ActionShortcut { 71 : public: 72 : explicit Angles(const ActionOptions&); 73 : static void registerKeywords( Keywords& keys ); 74 : }; 75 : 76 : PLUMED_REGISTER_ACTION(Angles,"ANGLES") 77 : 78 4 : void Angles::registerKeywords( Keywords& keys ) { 79 4 : ActionShortcut::registerKeywords( keys ); 80 4 : keys.add("atoms-1","GROUP","Calculate angles for each distinct set of three atoms in the group"); 81 4 : keys.add("atoms-2","GROUPA","A group of central atoms about which angles should be calculated"); 82 4 : keys.add("atoms-2","GROUPB","When used in conjunction with GROUPA this keyword instructs plumed " 83 : "to calculate all distinct angles involving one atom from GROUPA " 84 : "and two atoms from GROUPB. The atom from GROUPA is the central atom."); 85 4 : keys.add("atoms-3","GROUPC","This must be used in conjunction with GROUPA and GROUPB. All angles " 86 : "involving one atom from GROUPA, one atom from GROUPB and one atom from " 87 : "GROUPC are calculated. The GROUPA atoms are assumed to be the central " 88 : "atoms"); 89 4 : 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"); 90 4 : keys.addDOI("https://doi.org/10.1063/1.3628676"); 91 8 : keys.setValueDescription("vector","the ANGLE for each set of three atoms that were specified"); 92 4 : MultiColvarShortcuts::shortcutKeywords( keys ); 93 4 : keys.needsAction("ANGLE"); 94 4 : keys.needsAction("COORD_ANGLES"); 95 4 : } 96 : 97 1 : Angles::Angles(const ActionOptions&ao): 98 : Action(ao), 99 1 : ActionShortcut(ao) { 100 : std::string swit; 101 2 : parse("SWITCH",swit); 102 1 : if( swit.length()>0 ) { 103 : std::string cat, grp; 104 0 : parse("GROUPA",cat); 105 0 : parse("GROUPB",grp); 106 0 : if( cat.length()==0 || grp.length()==0 ) { 107 0 : error("must use GROUPA/GROUPB when using SWITCH"); 108 : } 109 0 : readInputLine( getShortcutLabel() + ": COORD_ANGLES SWITCH={" + swit + "} CATOMS=" + cat + " GROUP=" + grp + " " + convertInputLineToString() ); 110 : return; 111 : } 112 : std::vector<std::string> group; 113 2 : parseVector("GROUP",group); 114 : std::vector<std::string> groupa; 115 2 : parseVector("GROUPA",groupa); 116 : std::vector<std::string> groupb; 117 2 : parseVector("GROUPB",groupb); 118 : std::vector<std::string> groupc; 119 2 : parseVector("GROUPC",groupc); 120 1 : if( group.size()>0 ) { 121 0 : if( groupa.size()>0 || groupb.size()>0 || groupc.size()>0 ) { 122 0 : error("should only be GROUP keyword in input not GROUPA/GROUPB/GROUPC"); 123 : } 124 0 : Tools::interpretRanges( group ); 125 0 : std::string ainput = getShortcutLabel() + ": ANGLE"; 126 0 : unsigned n=1; 127 : // Not sure if this triple sum makes any sense 128 0 : for(unsigned i=2; i<group.size(); ++i ) { 129 0 : for(unsigned j=1; j<i; ++j ) { 130 0 : for(unsigned k=0; k<j; ++k) { 131 : std::string str_n; 132 0 : Tools::convert( n, str_n ); 133 0 : ainput += " ATOMS" + str_n + "=" + group[i] + "," + group[j] + "," + group[k]; 134 0 : n++; 135 : } 136 : } 137 : } 138 0 : readInputLine( ainput ); 139 1 : } else if( groupc.size()>0 ) { 140 0 : Tools::interpretRanges( groupa ); 141 0 : Tools::interpretRanges( groupb ); 142 0 : Tools::interpretRanges( groupc ); 143 0 : unsigned n=1; 144 0 : std::string ainput = getShortcutLabel() + ": ANGLE"; 145 0 : for(unsigned i=0; i<groupa.size(); ++i ) { 146 0 : for(unsigned j=0; j<groupb.size(); ++j ) { 147 0 : for(unsigned k=0; k<groupc.size(); ++k) { 148 : std::string str_n; 149 0 : Tools::convert( n, str_n ); 150 0 : ainput += " ATOMS" + str_n + "=" + groupb[j] + "," + groupa[i] + "," + groupc[k]; 151 0 : n++; 152 : } 153 : } 154 : } 155 0 : readInputLine( ainput ); 156 1 : } else if( groupa.size()>0 ) { 157 1 : Tools::interpretRanges( groupa ); 158 1 : Tools::interpretRanges( groupb ); 159 1 : unsigned n=1; 160 : std::string ainput; 161 1 : ainput = getShortcutLabel() + ": ANGLE"; 162 2 : for(unsigned i=0; i<groupa.size(); ++i ) { 163 9 : for(unsigned j=1; j<groupb.size(); ++j ) { 164 44 : for(unsigned k=0; k<j; ++k) { 165 : std::string str_n; 166 36 : Tools::convert( n, str_n ); 167 72 : ainput += " ATOMS" + str_n + "=" + groupb[j] + "," + groupa[i] + "," + groupb[k]; 168 36 : n++; 169 : } 170 : } 171 : } 172 1 : readInputLine( ainput ); 173 : } 174 1 : MultiColvarShortcuts::expandFunctions( getShortcutLabel(), getShortcutLabel(), "", this ); 175 1 : } 176 : 177 : } 178 : } 179 : 180 : 181 :