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 wham { 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. The following input demonstrates how this works in practise by showing you how this 34 : action can be used to analyze the output from a series of umbrella sampling calculations. 35 : The trajectory from each of the simulations run with the different biases should be concatenated into a 36 : single trajectory before running the following analysis script on the concatenated trajectory using PLUMED 37 : driver. The umbrella sampling simulations that will be analyzed using the script below applied a harmonic 38 : restraint that restrained the torsional angle involving atoms 5, 7, 9 and 15 to particular values. The script 39 : below calculates the reweighting weights for each of the trajectories and then applies the binless WHAM algorithm 40 : to determine a weight for each configuration in the concatenated trajectory. 41 : 42 : ```plumed 43 : #SETTINGS NREPLICAS=4 44 : phi: TORSION ATOMS=5,7,9,15 45 : rp: RESTRAINT ARG=phi KAPPA=50.0 AT=@replicas:{-3.00,-1.45,0.10,1.65} 46 : WHAM_WEIGHTS BIAS=rp.bias TEMP=300 FILE=wham-weights 47 : ``` 48 : 49 : The script above must be run with multiple replicas using the following command: 50 : 51 : ```` 52 : mpirun -np 4 plumed driver --mf_xtc alltraj.xtc --multi 4 53 : ```` 54 : 55 : For more details on how the weights for configurations are determined using the wham method see the documentation 56 : for the [WHAM](WHAM.md) action. 57 : 58 : */ 59 : //+ENDPLUMEDOC 60 : 61 : class WhamWeights : public ActionShortcut { 62 : public: 63 : static void registerKeywords( Keywords& keys ); 64 : explicit WhamWeights( const ActionOptions& ); 65 : }; 66 : 67 : PLUMED_REGISTER_ACTION(WhamWeights,"WHAM_WEIGHTS") 68 : 69 8 : void WhamWeights::registerKeywords( Keywords& keys ) { 70 8 : ActionShortcut::registerKeywords( keys ); 71 16 : keys.remove("LABEL"); 72 8 : keys.add("compulsory","BIAS","*.bias","the value of the biases to use when performing WHAM"); 73 8 : keys.add("optional","TEMP","the temperature at which the simulation was run"); 74 8 : keys.add("compulsory","STRIDE","1","the frequency with which the bias should be stored to perform WHAM"); 75 8 : keys.add("compulsory","FILE","the file on which to output the WHAM weights"); 76 8 : keys.add("optional","FMT","the format to use for the real numbers in the output file"); 77 16 : keys.setValueDescription("vector","the weights that were calculated using WHAM"); 78 8 : keys.needsAction("GATHER_REPLICAS"); 79 8 : keys.needsAction("CONCATENATE"); 80 8 : keys.needsAction("COLLECT"); 81 8 : keys.needsAction("WHAM"); 82 8 : keys.needsAction("DUMPVECTOR"); 83 8 : } 84 : 85 6 : WhamWeights::WhamWeights( const ActionOptions& ao ) : 86 : Action(ao), 87 6 : ActionShortcut(ao) { 88 : // Input for collection of weights for WHAM 89 : std::string bias; 90 12 : parse("BIAS",bias); 91 : std::string stride; 92 6 : parse("STRIDE",stride); 93 : // Input for GATHER_REPLICAS 94 12 : readInputLine( getShortcutLabel() + "_gather: GATHER_REPLICAS ARG=" + bias ); 95 : // Put all the replicas in a single vector 96 12 : readInputLine( getShortcutLabel() + "_gatherv: CONCATENATE ARG=" + getShortcutLabel() + "_gather.*"); 97 : // Input for COLLECT_FRAMES 98 12 : readInputLine( getShortcutLabel() + "_collect: COLLECT TYPE=vector ARG=" + getShortcutLabel() + "_gatherv STRIDE=" + stride); 99 : // Input for WHAM 100 6 : std::string temp, tempstr=""; 101 12 : parse("TEMP",temp); 102 6 : if( temp.length()>0 ) { 103 12 : tempstr="TEMP=" + temp; 104 : } 105 12 : readInputLine( getShortcutLabel() + ": WHAM ARG=" + getShortcutLabel() + "_collect " + tempstr ); 106 : // Input for PRINT (will just output at end of calc 107 : std::string filename, fmt; 108 6 : parse("FILE",filename); 109 12 : parse("FMT",fmt); 110 6 : if(fmt.length()>0) { 111 12 : fmt=" FMT=" + fmt; 112 : } 113 12 : readInputLine( "DUMPVECTOR STRIDE=0 ARG=" + getShortcutLabel() + " FILE=" + filename + fmt ); 114 6 : } 115 : 116 : } 117 : }