Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2017-2020 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/ActionRegister.h" 23 : #include "core/ActionWithValue.h" 24 : #include "core/ActionWithArguments.h" 25 : #include "core/PlumedMain.h" 26 : #include "tools/Communicator.h" 27 : 28 : //+PLUMEDOC REWEIGHTING WHAM 29 : /* 30 : Calculate the weights for configurations using the weighted histogram analysis method. 31 : 32 : Suppose that you have run multiple \f$N\f$ trajectories each of which was computed by integrating a different biased Hamiltonian. We can calculate the probability of observing 33 : the set of configurations during the \f$N\f$ trajectories that we ran using the product of multinomial distributions shown below: 34 : \f[ 35 : P( \vec{T} ) \propto \prod_{j=1}^M \prod_{k=1}^N (c_k w_{kj} p_j)^{t_{kj}} 36 : \label{eqn:wham1} 37 : \f] 38 : In this expression the second product runs over the biases that were used when calculating the \f$N\f$ trajectories. The first product then runs over the 39 : \f$M\f$ bins in our histogram. The \f$p_j\f$ variable that is inside this product is the quantity we wish to extract; namely, the unbiased probability of 40 : having a set of CV values that lie within the range for the \f$j\f$th bin. 41 : 42 : The quantity that we can easily extract from our simulations, \f$t_{kj}\f$, however, measures the number of frames from trajectory \f$k\f$ that are inside the \f$j\f$th bin. 43 : To interpret this quantity we must consider the bias that acts on each of the replicas so the \f$w_{kj}\f$ term is introduced. This quantity is calculated as: 44 : \f[ 45 : w_{kj} = 46 : \f] 47 : and is essentially the factor that we have to multiply the unbiased probability of being in the bin by in order to get the probability that we would be inside this same bin in 48 : the \f$k\f$th of our biased simulations. Obviously, these \f$w_{kj}\f$ values depend on the value that the CVs take and also on the particular trajectory that we are investigating 49 : all of which, remember, have different simulation biases. Finally, \f$c_k\f$ is a free parameter that ensures that, for each \f$k\f$, the biased probability is normalized. 50 : 51 : We can use the equation for the probability that was given above to find a set of values for \f$p_j\f$ that maximizes the likelihood of observing the trajectory. 52 : This constrained optimization must be performed using a set of Lagrange multipliers, \f$\lambda_k\f$, that ensure that each of the biased probability distributions 53 : are normalized, that is \f$\sum_j c_kw_{kj}p_j=1\f$. Furthermore, as the problem is made easier if the quantity in the equation above is replaced by its logarithm 54 : we actually chose to minimize 55 : \f[ 56 : \mathcal{L}= \sum_{j=1}^M \sum_{k=1}^N t_{kj} \ln c_k w_{kj} p_j + \sum_k\lambda_k \left( \sum_{j=1}^N c_k w_{kj} p_j - 1 \right) 57 : \f] 58 : After some manipulations the following (WHAM) equations emerge: 59 : \f[ 60 : \begin{aligned} 61 : p_j & \propto \frac{\sum_{k=1}^N t_{kj}}{\sum_k c_k w_{kj}} \\ 62 : c_k & =\frac{1}{\sum_{j=1}^M w_{kj} p_j} 63 : \end{aligned} 64 : \f] 65 : which can be solved by computing the \f$p_j\f$ values using the first of the two equations above with an initial guess for the \f$c_k\f$ values and by then refining 66 : these \f$p_j\f$ values using the \f$c_k\f$ values that are obtained by inserting the \f$p_j\f$ values obtained into the second of the two equations above. 67 : 68 : Notice that only \f$\sum_k t_{kj}\f$, which is the total number of configurations from all the replicas that enter the \f$j\f$th bin, enters the WHAM equations above. 69 : There is thus no need to record which replica generated each of the frames. One can thus simply gather the trajectories from all the replicas together at the outset. 70 : This observation is important as it is the basis of the binless formulation of WHAM that is implemented within PLUMED. 71 : 72 : \par Examples 73 : 74 : */ 75 : //+ENDPLUMEDOC 76 : 77 : namespace PLMD { 78 : namespace wham { 79 : 80 : class Wham : 81 : public ActionWithValue, 82 : public ActionWithArguments { 83 : private: 84 : double thresh, simtemp; 85 : unsigned nreplicas; 86 : unsigned maxiter; 87 : public: 88 : static void registerKeywords(Keywords&); 89 : explicit Wham(const ActionOptions&ao); 90 0 : unsigned getNumberOfDerivatives() { return 0; } 91 : void calculate() override ; 92 0 : void apply() override {} 93 : }; 94 : 95 : PLUMED_REGISTER_ACTION(Wham,"WHAM") 96 : 97 29 : void Wham::registerKeywords(Keywords& keys ) { 98 29 : Action::registerKeywords( keys ); ActionWithValue::registerKeywords( keys ); 99 29 : ActionWithArguments::registerKeywords( keys ); 100 58 : keys.addInputKeyword("compulsory","ARG","scalar/vector/matrix","the stored values for the bias"); 101 58 : keys.add("compulsory","MAXITER","1000","maximum number of iterations for WHAM algorithm"); 102 58 : keys.add("compulsory","WHAMTOL","1e-10","threshold for convergence of WHAM algorithm"); 103 58 : keys.add("optional","TEMP","the system temperature. This is not required if your MD code passes this quantity to PLUMED"); 104 29 : keys.remove("NUMERICAL_DERIVATIVES"); 105 58 : keys.setValueDescription("vector","the vector of WHAM weights to use for reweighting the elements in a time series"); 106 29 : } 107 : 108 13 : Wham::Wham(const ActionOptions&ao): 109 : Action(ao), 110 : ActionWithValue(ao), 111 13 : ActionWithArguments(ao) 112 : { 113 : // Read in the temperature 114 13 : simtemp=getkBT(); 115 13 : if(simtemp==0) error("The MD engine does not pass the temperature to plumed so you have to specify it using TEMP"); 116 : // Now read in parameters of WHAM 117 26 : parse("MAXITER",maxiter); parse("WHAMTOL",thresh); 118 13 : if(comm.Get_rank()==0) nreplicas=multi_sim_comm.Get_size(); 119 13 : comm.Bcast(nreplicas,0); addValue( getPntrToArgument(0)->getShape() ); setNotPeriodic(); 120 13 : } 121 : 122 12 : void Wham::calculate() { 123 : // Retrieve the values that were stored for the biase 124 12 : std::vector<double> stored_biases( getPntrToArgument(0)->getNumberOfValues() ); 125 25284 : for(unsigned i=0; i<stored_biases.size(); ++i) stored_biases[i] = getPntrToArgument(0)->get(i); 126 : // Get the minimum value of the bias 127 12 : double minv = *min_element(std::begin(stored_biases), std::end(stored_biases)); 128 : // Resize final weights array 129 12 : plumed_assert( stored_biases.size()%nreplicas==0 ); 130 12 : std::vector<double> final_weights( stored_biases.size() / nreplicas, 1.0 ); 131 12 : if( getPntrToComponent(0)->getNumberOfValues()!=final_weights.size() ) { 132 12 : std::vector<unsigned> shape(1); shape[0]=final_weights.size(); getPntrToComponent(0)->setShape( shape ); 133 : } 134 : // Offset and exponential of the bias 135 12 : std::vector<double> expv( stored_biases.size() ); 136 25284 : for(unsigned i=0; i<expv.size(); ++i) expv[i] = exp( (-stored_biases[i]+minv) / simtemp ); 137 : // Initialize Z 138 12 : std::vector<double> Z( nreplicas, 1.0 ), oldZ( nreplicas ); 139 : // Now the iterative loop to calculate the WHAM weights 140 4584 : for(unsigned iter=0; iter<maxiter; ++iter) { 141 : // Store Z 142 32088 : for(unsigned j=0; j<Z.size(); ++j) oldZ[j]=Z[j]; 143 : // Recompute weights 144 : double norm=0; 145 1613568 : for(unsigned j=0; j<final_weights.size(); ++j) { 146 : double ew=0; 147 11262888 : for(unsigned k=0; k<Z.size(); ++k) ew += expv[j*Z.size()+k] / Z[k]; 148 1608984 : final_weights[j] = 1.0 / ew; norm += final_weights[j]; 149 : } 150 : // Normalize weights 151 1613568 : for(unsigned j=0; j<final_weights.size(); ++j) final_weights[j] /= norm; 152 : // Recompute Z 153 32088 : for(unsigned j=0; j<Z.size(); ++j) Z[j] = 0.0; 154 1613568 : for(unsigned j=0; j<final_weights.size(); ++j) { 155 11262888 : for(unsigned k=0; k<Z.size(); ++k) Z[k] += final_weights[j]*expv[j*Z.size()+k]; 156 : } 157 : // Normalize Z and compute change in Z 158 32088 : double change=0; norm=0; for(unsigned k=0; k<Z.size(); ++k) norm+=Z[k]; 159 32088 : for(unsigned k=0; k<Z.size(); ++k) { 160 27504 : Z[k] /= norm; double d = std::log( Z[k] / oldZ[k] ); change += d*d; 161 : } 162 4584 : if( change<thresh ) { 163 4224 : for(unsigned j=0; j<final_weights.size(); ++j) getPntrToComponent(0)->set( j, final_weights[j] ); 164 12 : return; 165 : } 166 : } 167 0 : error("Too many iterations in WHAM" ); 168 : } 169 : 170 : } 171 : }