Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2011-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 "Colvar.h" 23 : #include "ActionRegister.h" 24 : #include "core/PlumedMain.h" 25 : #include "core/Atoms.h" 26 : 27 : #include <string> 28 : #include <cmath> 29 : 30 : namespace PLMD { 31 : namespace colvar { 32 : 33 : //+PLUMEDOC COLVAR ENERGY 34 : /* 35 : Calculate the total potential energy of the simulation box. 36 : 37 : The potential energy can be biased e.g. with umbrella sampling \cite bart-karp98jpcb or with well-tempered metadynamics \cite Bonomi:2009p17935. 38 : 39 : Notice that this CV could be unavailable with some MD code. When 40 : it is available, and when also replica exchange is available, 41 : metadynamics applied to ENERGY can be used to decrease the 42 : number of required replicas. 43 : 44 : \bug This \ref ENERGY does not include long tail corrections. 45 : Thus when using e.g. LAMMPS `"pair_modify tail yes"` or GROMACS `"DispCorr Ener"` (or `"DispCorr EnerPres"`), 46 : the potential energy from \ref ENERGY will be slightly different form the one of the MD code. 47 : You should still be able to use \ref ENERGY and then reweight your simulation with the correct MD energy value. 48 : 49 : \bug Acceptance for replica exchange when \ref ENERGY is biased 50 : is computed correctly only if all the replicas have the same 51 : potential energy function. This is for instance not true when 52 : using GROMACS with lambda replica exchange or with plumed-hrex branch. 53 : 54 : \par Examples 55 : 56 : The following input instructs plumed to print the energy of the system 57 : \plumedfile 58 : ene: ENERGY 59 : PRINT ARG=ene 60 : \endplumedfile 61 : 62 : */ 63 : //+ENDPLUMEDOC 64 : 65 : 66 78 : class Energy : public Colvar { 67 : 68 : public: 69 : explicit Energy(const ActionOptions&); 70 : // active methods: 71 : void prepare(); 72 : virtual void calculate(); 73 : unsigned getNumberOfDerivatives(); 74 : static void registerKeywords( Keywords& keys ); 75 : }; 76 : 77 : 78 : using namespace std; 79 : 80 : 81 7434 : PLUMED_REGISTER_ACTION(Energy,"ENERGY") 82 : 83 39 : Energy::Energy(const ActionOptions&ao): 84 39 : PLUMED_COLVAR_INIT(ao) 85 : { 86 : // if(checkNumericalDerivatives()) 87 : // error("Cannot use NUMERICAL_DERIVATIVES with ENERGY"); 88 39 : isEnergy=true; 89 39 : addValueWithDerivatives(); setNotPeriodic(); 90 39 : getPntrToValue()->resizeDerivatives(1); 91 39 : log<<" Bibliography "; 92 117 : log<<plumed.cite("Bartels and Karplus, J. Phys. Chem. B 102, 865 (1998)"); 93 117 : log<<plumed.cite("Bonomi and Parrinello, J. Comp. Chem. 30, 1615 (2009)"); 94 39 : log<<"\n"; 95 39 : } 96 : 97 40 : void Energy::registerKeywords( Keywords& keys ) { 98 40 : Action::registerKeywords( keys ); 99 40 : ActionAtomistic::registerKeywords( keys ); 100 40 : ActionWithValue::registerKeywords( keys ); 101 80 : keys.remove("NUMERICAL_DERIVATIVES"); 102 40 : } 103 : 104 2 : unsigned Energy::getNumberOfDerivatives() { 105 2 : return 1; 106 : } 107 : 108 3888 : void Energy::prepare() { 109 3888 : plumed.getAtoms().setCollectEnergy(true); 110 3888 : } 111 : 112 : // calculator 113 3888 : void Energy::calculate() { 114 7776 : setValue( getEnergy() ); 115 3888 : getPntrToComponent(0)->addDerivative(0,1.0); 116 3888 : } 117 : 118 : } 119 5517 : } 120 : 121 : 122 :