All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Min.cpp
Go to the documentation of this file.
1 /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2  Copyright (c) 2013 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-code.org for more information.
6 
7  This file is part of plumed, version 2.0.
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 "VesselRegister.h"
23 #include "FunctionVessel.h"
24 #include "ActionWithVessel.h"
25 
26 namespace PLMD {
27 namespace vesselbase{
28 
29 class Min : public FunctionVessel {
30 private:
31  std::vector<double> df;
32  double beta;
33 public:
34  static void registerKeywords( Keywords& keys );
35  static void reserveKeyword( Keywords& keys );
36  Min( const VesselOptions& da );
37  std::string function_description();
38  bool calculate();
39  void finish();
40 };
41 
43 
44 void Min::registerKeywords( Keywords& keys ){
46  keys.add("compulsory","BETA","the value of beta for the equation in the manual");
47 }
48 
50  keys.reserve("optional","MIN","calculate the minimum value. "
51  "To make this quantity continuous the minimum is calculated using "
52  "\\f$ \\textrm{min} = \\frac{\\beta}{ \\log \\sum_i \\exp\\left( \\frac{\\beta}{s_i} \\right) } \\f$ "
53  "The value of \\f$\\beta\\f$ in this function is specified using (BETA=\\f$\\beta\\f$)",true);
54  keys.addOutputComponent("min","MIN","the minimum value. This is calculated using the formula described in the description of the "
55  "keyword so as to make it continuous.");
56 }
57 
59 FunctionVessel(da),
60 df(2)
61 {
62  if( getAction()->isPeriodic() ) error("min is not a meaningful option for periodic variables");
63  parse("BETA",beta);
64 
65  if( diffweight ) error("can't calculate min if weight is differentiable");
66 }
67 
69  std::string str_beta; Tools::convert( beta, str_beta );
70  return "the minimum value. Beta is equal to " + str_beta;
71 }
72 
74  double val=getAction()->getElementValue(0);
75  double dval, f = exp(beta/val); dval=f/(val*val);
77  getAction()->chainRuleForElementDerivatives( 0, 0, dval, this );
78  return true;
79 }
80 
81 void Min::finish(){
82  double valin=getFinalValue(0); double dist=beta/std::log( valin );
83  setOutputValue( dist ); df[0]=dist*dist/valin; df[1]=0.0;
85 }
86 
87 }
88 }
double getFinalValue(const unsigned &j)
Get the nth value in the distribution.
static void reserveKeyword(Keywords &keys)
Definition: Min.cpp:49
bool diffweight
Are the derivatives differentiable.
void chainRuleForElementDerivatives(const unsigned &, const unsigned &, const double &, Vessel *)
Merge the derivatives.
static bool convert(const std::string &str, double &t)
Convert a string to a double, reading it.
Definition: Tools.cpp:74
void setOutputValue(const double &val)
Set the final value.
Objects that inherit from FunctionVessel can be used (in tandem with PLMD::ActionWithVessel) to calcu...
static void registerKeywords(Keywords &keys)
void mergeFinalDerivatives(const std::vector< double > &df)
This does a combination of the product and chain rules.
bool calculate()
Calculate the part of the vessel that is done in the loop.
Definition: Min.cpp:73
void addOutputComponent(const std::string &name, const std::string &key, const std::string &descr)
Add a potential component which can be output by this particular action.
Definition: Keywords.cpp:565
This class holds the keywords and their documentation.
Definition: Keywords.h:36
void addValueIgnoringTolerance(const unsigned &jval, const double &val)
Add some value to the accumulator and ignore the tolerance.
#define PLUMED_REGISTER_VESSEL(classname, keyword)
Min(const VesselOptions &da)
Definition: Min.cpp:58
void reserve(const std::string &t, const std::string &k, const std::string &d, const bool isvessel=false)
Reserve a keyword.
Definition: Keywords.cpp:110
void error(const std::string &errmsg)
Report an error.
Definition: Vessel.cpp:122
This class is used to pass the input to Vessels.
Definition: Vessel.h:53
std::vector< double > df
Definition: Min.cpp:31
void int double * da
Definition: Matrix.h:47
void parse(const std::string &key, T &t)
Parse something from the input.
Definition: Vessel.h:169
void finish()
Complete the calculation once the loop is finished.
Definition: Min.cpp:81
static void registerKeywords(Keywords &keys)
Definition: Min.cpp:44
std::string function_description()
The rest of the description of what we are calculating.
Definition: Min.cpp:68
double getElementValue(const unsigned &ival) const
Get the value of this element.
ActionWithVessel * getAction()
Return a pointer to the action we are working in.
Definition: Vessel.h:236