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/ActionRegister.h" 24 : 25 : namespace PLMD { 26 : namespace analysis { 27 : 28 : //+PLUMEDOC REWEIGHTING WHAM_WEIGHTS 29 : /* 30 : Calculate and output weights for configurations using the weighted histogram analysis method. 31 : 32 : This shortcut action allows you to calculate and output weights computed using the weighted histogram 33 : analysis technique. For more detail on how this technique works see \ref REWEIGHT_WHAM 34 : 35 : \par Examples 36 : 37 : The following input can be used to analyze the output from a series of umbrella sampling calculations. 38 : The trajectory from each of the simulations run with the different biases should be concatenated into a 39 : single trajectory before running the following analysis script on the concatenated trajectory using PLUMED 40 : driver. The umbrella sampling simulations that will be analyzed using the script below applied a harmonic 41 : restraint that restrained the torsional angle involving atoms 5, 7, 9 and 15 to particular values. The script 42 : below calculates the reweighting weights for each of the trajectories and then applies the binless WHAM algorithm 43 : to determine a weight for each configuration in the concatenated trajectory. 44 : 45 : \plumedfile 46 : #SETTINGS NREPLICAS=4 47 : phi: TORSION ATOMS=5,7,9,15 48 : rp: RESTRAINT ARG=phi KAPPA=50.0 ... 49 : AT=@replicas:{ 50 : -3.00000000000000000000 51 : -1.45161290322580645168 52 : .09677419354838709664 53 : 1.64516129032258064496 54 : } 55 : ... 56 : 57 : WHAM_WEIGHTS BIAS=rp.bias TEMP=300 FILE=wham-weights 58 : \endplumedfile 59 : 60 : The script above must be run with multiple replicas using the following command: 61 : 62 : \verbatim 63 : mpirun -np 4 plumed driver --mf_xtc alltraj.xtc --multi 4 64 : \endverbatim 65 : 66 : */ 67 : //+ENDPLUMEDOC 68 : 69 : class WhamWeights : public ActionShortcut { 70 : public: 71 : static void registerKeywords( Keywords& keys ); 72 : explicit WhamWeights( const ActionOptions& ); 73 : }; 74 : 75 10431 : PLUMED_REGISTER_ACTION(WhamWeights,"WHAM_WEIGHTS") 76 : 77 7 : void WhamWeights::registerKeywords( Keywords& keys ) { 78 7 : ActionShortcut::registerKeywords( keys ); keys.remove("LABEL"); 79 14 : keys.add("compulsory","BIAS","*.bias","the value of the biases to use when performing WHAM"); 80 14 : keys.add("compulsory","TEMP","the temperature at which the simulation was run"); 81 14 : keys.add("compulsory","STRIDE","1","the frequency with which the bias should be stored to perform WHAM"); 82 14 : keys.add("compulsory","FILE","the file on which to output the WHAM weights"); 83 14 : keys.add("optional","FMT","the format to use for the real numbers in the output file"); 84 7 : } 85 : 86 6 : WhamWeights::WhamWeights( const ActionOptions& ao ) : 87 : Action(ao), 88 6 : ActionShortcut(ao) 89 : { 90 : // Input for REWEIGHT_WHAM 91 6 : std::string rew_line = getShortcutLabel() + "_weights: REWEIGHT_WHAM"; 92 18 : std::string bias; parse("BIAS",bias); rew_line += " ARG=" + bias; 93 12 : std::string temp; parse("TEMP",temp); rew_line += " TEMP=" + temp; 94 6 : readInputLine( rew_line ); 95 : // Input for COLLECT_FRAMES 96 12 : std::string col_line = getShortcutLabel() + "_collect: COLLECT_FRAMES LOGWEIGHTS=" + getShortcutLabel() + "_weights"; 97 12 : std::string stride; parse("STRIDE",stride); col_line += " STRIDE=" + stride; 98 6 : readInputLine( col_line ); 99 : // Input for line to output data 100 12 : std::string out_line="OUTPUT_ANALYSIS_DATA_TO_COLVAR USE_OUTPUT_DATA_FROM=" + getShortcutLabel() + "_collect"; 101 12 : std::string file; parse("FILE",file); out_line += " FILE=" + file; 102 12 : std::string fmt="%f"; parse("FMT",fmt); out_line += " FMT=" + fmt; 103 6 : readInputLine( out_line ); 104 6 : } 105 : 106 : } 107 : }