Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2011-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/ActionToPutData.h" 23 : #include "core/DomainDecomposition.h" 24 : #include "core/PlumedMain.h" 25 : #include "core/ActionSet.h" 26 : #include "core/ActionRegister.h" 27 : 28 : namespace PLMD { 29 : namespace colvar { 30 : 31 : //+PLUMEDOC COLVAR ENERGY 32 : /* 33 : Calculate the total potential energy of the simulation box. 34 : 35 : As is explained in the papers in the bibliography the potential energy can be biased with umbrella sampling. 36 : To print the potential energy from PLUMED you can use an input similar to the one below: 37 : 38 : ```plumed 39 : ene: ENERGY 40 : PRINT ARG=ene FILE=colvar 41 : ``` 42 : 43 : Notice that this CV is not available with all the MD codes. When 44 : it is available, and when replica exchange is also available, 45 : metadynamics applied to ENERGY can be used to decrease the 46 : number of required replicas. 47 : 48 : > [!CAUTION] 49 : > The ENERGY output by PLUMED does not include long tail corrections. 50 : > Thus when using e.g. LAMMPS `"pair_modify tail yes"` or GROMACS `"DispCorr Ener"` (or `"DispCorr EnerPres"`), 51 : > the potential energy from ENERGY will be slightly different from the one that is output by the MD code. 52 : > You should still be able to bias the ENERGY and then reweight your simulation with the correct MD energy values. 53 : 54 : > [!CAUTION] 55 : > Acceptance for replica exchange when ENERGY is biased 56 : > is computed correctly only if all the replicas have the same 57 : > potential energy function. This is for instance not true when 58 : > using GROMACS with lambda replica exchange or with plumed-hrex branch. 59 : 60 : */ 61 : //+ENDPLUMEDOC 62 : 63 : 64 : class Energy : public ActionToPutData { 65 : private: 66 : /// This is used to sum the data 67 : DomainDecomposition* interface; 68 : /// This is the list of forces that must be scaled 69 : std::vector<ActionToPutData*> forces_to_scale; 70 : public: 71 : explicit Energy(const ActionOptions&); 72 : // active methods: 73 : static void registerKeywords( Keywords& keys ); 74 : void wait() override; 75 : void apply() override; 76 : }; 77 : 78 : 79 : PLUMED_REGISTER_ACTION(Energy,"ENERGY") 80 : 81 40 : Energy::Energy(const ActionOptions&ao): 82 : Action(ao), 83 : ActionToPutData(ao), 84 40 : interface(NULL) { 85 40 : plumed.setEnergyValue( getLabel() ); 86 : std::vector<unsigned> shape; 87 40 : addValue( shape ); 88 40 : setNotPeriodic(); 89 80 : setUnit( "energy", "default" ); 90 40 : ActionToPutData* px=plumed.getActionSet().selectWithLabel< ActionToPutData*>("posx"); 91 40 : plumed_assert(px); 92 40 : forces_to_scale.push_back(px); 93 40 : addDependency( px ); 94 40 : ActionToPutData* py=plumed.getActionSet().selectWithLabel< ActionToPutData*>("posy"); 95 40 : plumed_assert(py); 96 40 : forces_to_scale.push_back(py); 97 40 : addDependency( py ); 98 40 : ActionToPutData* pz=plumed.getActionSet().selectWithLabel< ActionToPutData*>("posz"); 99 40 : plumed_assert(pz); 100 40 : forces_to_scale.push_back(pz); 101 40 : addDependency( pz ); 102 40 : ActionToPutData* bx=plumed.getActionSet().selectWithLabel< ActionToPutData*>("Box"); 103 40 : plumed_assert(bx); 104 40 : forces_to_scale.push_back(bx); 105 40 : addDependency( bx ); 106 40 : log<<" Bibliography "; 107 80 : log<<plumed.cite("Bartels and Karplus, J. Phys. Chem. B 102, 865 (1998)"); 108 80 : log<<plumed.cite("Bonomi and Parrinello, J. Comp. Chem. 30, 1615 (2009)"); 109 40 : log<<"\n"; 110 40 : } 111 : 112 42 : void Energy::registerKeywords( Keywords& keys ) { 113 42 : Action::registerKeywords( keys ); 114 84 : keys.setValueDescription("scalar","the internal energy"); 115 42 : keys.addDOI("10.1021/jp972280j"); 116 42 : keys.addDOI("10.1103/PhysRevLett.104.190601"); 117 42 : } 118 : 119 3989 : void Energy::wait() { 120 3989 : if( !interface ) { 121 40 : std::vector<DomainDecomposition*> allput=plumed.getActionSet().select<DomainDecomposition*>(); 122 40 : if( allput.size()>1 ) { 123 0 : warning("found more than one interface so don't know how to sum energy"); 124 : } 125 40 : interface = allput[0]; 126 : } 127 3989 : ActionToPutData::wait(); 128 3989 : if( interface ) { 129 3989 : interface->sumOverDomains( copyOutput(0) ); 130 : } 131 3989 : } 132 : 133 3989 : void Energy::apply() { 134 3989 : if( getPntrToValue()->forcesWereAdded() ) { 135 250 : for(unsigned i=0; i<forces_to_scale.size(); ++i) { 136 200 : forces_to_scale[i]->rescaleForces( 1.- getPntrToValue()->getForce(0)); 137 : } 138 : } 139 3989 : } 140 : 141 : } 142 : } 143 : 144 : 145 :