Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2018-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/PlumedMain.h" 24 : #include "core/ActionSet.h" 25 : #include "core/ActionRegister.h" 26 : 27 : namespace PLMD { 28 : namespace wham { 29 : 30 : //+PLUMEDOC REWEIGHTING WHAM_HISTOGRAM 31 : /* 32 : This can be used to output the a histogram using the weighted histogram technique 33 : 34 : This shortcut action allows you to calculate a histogram using the weighted histogram analysis technique. 35 : The following input illustrates how this is used in practise to analyze the output from a series of umbrella sampling calculations. 36 : The trajectory from each of the simulations run with the different biases should be concatenated into a 37 : single trajectory before running the following analysis script on the concatenated trajectory using PLUMED 38 : driver. The umbrella sampling simulations that will be analyzed using the script below applied a harmonic 39 : restraint that restrained the torsional angle involving atoms 5, 7, 9 and 15 to particular values. The script 40 : below calculates the reweighting weights for each of the trajectories and then applies the binless WHAM algorithm 41 : to determine a weight for each configuration in the concatenated trajectory. A histogram is then constructed from 42 : the configurations visited and their weights. This histogram is then converted into a free energy surface and output 43 : to a file called fes.dat 44 : 45 : ```plumed 46 : #SETTINGS NREPLICAS=4 47 : phi: TORSION ATOMS=5,7,9,15 48 : psi: TORSION ATOMS=7,9,15,17 49 : rp: RESTRAINT ARG=phi KAPPA=50.0 AT=@replicas:{-3.00,-1.45,0.10,1.65} 50 : hh: WHAM_HISTOGRAM ARG=phi BIAS=rp.bias TEMP=300 GRID_MIN=-pi GRID_MAX=pi GRID_BIN=50 51 : fes: CONVERT_TO_FES ARG=hh TEMP=300 52 : DUMPGRID ARG=fes FILE=fes.dat 53 : ``` 54 : 55 : The script above must be run with multiple replicas using the following command: 56 : 57 : ```` 58 : mpirun -np 4 plumed driver --mf_xtc alltraj.xtc --multi 4 59 : ```` 60 : 61 : For more details on how the weights for configurations are determined using the wham method see the documentation 62 : for the [WHAM](WHAM.md) action. 63 : 64 : */ 65 : //+ENDPLUMEDOC 66 : 67 : class WhamHistogram : public ActionShortcut { 68 : public: 69 : static void registerKeywords( Keywords& keys ); 70 : explicit WhamHistogram( const ActionOptions& ); 71 : }; 72 : 73 : PLUMED_REGISTER_ACTION(WhamHistogram,"WHAM_HISTOGRAM") 74 : 75 10 : void WhamHistogram::registerKeywords( Keywords& keys ) { 76 10 : ActionShortcut::registerKeywords( keys ); 77 10 : keys.add("compulsory","ARG","the arguments that you would like to make the histogram for"); 78 10 : keys.add("compulsory","BIAS","*.bias","the value of the biases to use when performing WHAM"); 79 10 : keys.add("compulsory","TEMP","the temperature at which the simulation was run"); 80 10 : keys.add("compulsory","STRIDE","1","the frequency with which the data should be stored to perform WHAM"); 81 10 : keys.add("compulsory","GRID_MIN","the minimum to use for the grid"); 82 10 : keys.add("compulsory","GRID_MAX","the maximum to use for the grid"); 83 10 : keys.add("compulsory","GRID_BIN","the number of bins to use for the grid"); 84 10 : keys.add("optional","BANDWIDTH","the bandwidth for kernel density estimation"); 85 20 : keys.setValueDescription("grid","the histogram that was generated using the WHAM weights"); 86 10 : keys.needsAction("GATHER_REPLICAS"); 87 10 : keys.needsAction("CONCATENATE"); 88 10 : keys.needsAction("COLLECT"); 89 10 : keys.needsAction("WHAM"); 90 10 : keys.needsAction("KDE"); 91 10 : } 92 : 93 : 94 7 : WhamHistogram::WhamHistogram( const ActionOptions& ao ) : 95 : Action(ao), 96 7 : ActionShortcut(ao) { 97 : // Input for collection of weights for WHAM 98 : std::string bias; 99 14 : parse("BIAS",bias); 100 : std::string stride; 101 7 : parse("STRIDE",stride); 102 : // Input for GATHER_REPLICAS 103 14 : readInputLine( getShortcutLabel() + "_gather: GATHER_REPLICAS ARG=" + bias ); 104 : // Put all the replicas in a single vector 105 14 : readInputLine( getShortcutLabel() + "_gatherv: CONCATENATE ARG=" + getShortcutLabel() + "_gather.*"); 106 : // Input for COLLECT_FRAMES 107 14 : readInputLine( getShortcutLabel() + "_collect: COLLECT TYPE=vector ARG=" + getShortcutLabel() + "_gatherv STRIDE=" + stride); 108 : // Input for WHAM 109 7 : std::string temp, tempstr=""; 110 14 : parse("TEMP",temp); 111 7 : if( temp.length()>0 ) { 112 14 : tempstr="TEMP=" + temp; 113 : } 114 14 : readInputLine( getShortcutLabel() + "_wham: WHAM ARG=" + getShortcutLabel() + "_collect " + tempstr ); 115 : // Input for COLLECT_FRAMES 116 : std::vector<std::string> args; 117 14 : parseVector("ARG",args); 118 : std::string argstr; 119 14 : for(unsigned i=0; i<args.size(); ++i) { 120 14 : readInputLine( getShortcutLabel() + "_data_" + args[i] + ": COLLECT ARG=" + args[i] ); 121 7 : if( i==0 ) { 122 : argstr = " ARG="; 123 : } else { 124 : argstr += ","; 125 : } 126 14 : argstr += getShortcutLabel() + "_data_" + args[i]; 127 : } 128 : // Input for HISTOGRAM 129 7 : std::string histo_line, bw=""; 130 14 : parse("BANDWIDTH",bw); 131 7 : if( bw!="" ) { 132 0 : histo_line += " BANDWIDTH=" + bw; 133 : } else { 134 : histo_line += " KERNEL=DISCRETE"; 135 : } 136 : std::string min; 137 7 : parse("GRID_MIN",min); 138 14 : histo_line += " GRID_MIN=" + min; 139 : std::string max; 140 7 : parse("GRID_MAX",max); 141 14 : histo_line += " GRID_MAX=" + max; 142 : std::string bin; 143 7 : parse("GRID_BIN",bin); 144 7 : histo_line += " GRID_BIN=" + bin; 145 14 : readInputLine( getShortcutLabel() + ": KDE " + argstr + " HEIGHTS=" + getShortcutLabel() + "_wham" + histo_line ); 146 14 : } 147 : 148 : } 149 : }