Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2011-2019 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 energy of the simulation box. 36 : 37 : Total energy can be biased 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 Acceptance for replica exchange when \ref ENERGY is biased 45 : is computed correctly only of all the replicas has the same 46 : potential energy function. This is for instance not true when 47 : using GROMACS with lambda replica exchange of with plumed-hrex branch. 48 : 49 : \par Examples 50 : 51 : The following input instructs plumed to print the energy of the system 52 : \plumedfile 53 : ene: ENERGY 54 : PRINT ARG=ene 55 : \endplumedfile 56 : 57 : */ 58 : //+ENDPLUMEDOC 59 : 60 : 61 76 : class Energy : public Colvar { 62 : 63 : public: 64 : explicit Energy(const ActionOptions&); 65 : // active methods: 66 : void prepare() override; 67 : void calculate() override; 68 : unsigned getNumberOfDerivatives() override; 69 : static void registerKeywords( Keywords& keys ); 70 : }; 71 : 72 : 73 : using namespace std; 74 : 75 : 76 7908 : PLUMED_REGISTER_ACTION(Energy,"ENERGY") 77 : 78 38 : Energy::Energy(const ActionOptions&ao): 79 38 : PLUMED_COLVAR_INIT(ao) 80 : { 81 : // if(checkNumericalDerivatives()) 82 : // error("Cannot use NUMERICAL_DERIVATIVES with ENERGY"); 83 38 : isEnergy=true; 84 38 : addValueWithDerivatives(); setNotPeriodic(); 85 38 : getPntrToValue()->resizeDerivatives(1); 86 38 : log<<" Bibliography "; 87 114 : log<<plumed.cite("Bartels and Karplus, J. Phys. Chem. B 102, 865 (1998)"); 88 114 : log<<plumed.cite("Bonomi and Parrinello, J. Comp. Chem. 30, 1615 (2009)"); 89 38 : log<<"\n"; 90 38 : } 91 : 92 39 : void Energy::registerKeywords( Keywords& keys ) { 93 39 : Action::registerKeywords( keys ); 94 39 : ActionAtomistic::registerKeywords( keys ); 95 39 : ActionWithValue::registerKeywords( keys ); 96 78 : keys.remove("NUMERICAL_DERIVATIVES"); 97 39 : } 98 : 99 2 : unsigned Energy::getNumberOfDerivatives() { 100 2 : return 1; 101 : } 102 : 103 3787 : void Energy::prepare() { 104 3787 : plumed.getAtoms().setCollectEnergy(true); 105 3787 : } 106 : 107 : // calculator 108 3787 : void Energy::calculate() { 109 3787 : setValue( getEnergy() ); 110 3787 : getPntrToComponent(0)->addDerivative(0,1.0); 111 3787 : } 112 : 113 : } 114 5874 : } 115 : 116 : 117 :