Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2016-2018 The VES code team 3 : (see the PEOPLE-VES file at the root of this folder for a list of names) 4 : 5 : See http://www.ves-code.org for more information. 6 : 7 : This file is part of VES code module. 8 : 9 : The VES code module 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 : The VES code module 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 the VES code module. If not, see <http://www.gnu.org/licenses/>. 21 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */ 22 : 23 : #include "CoeffsVector.h" 24 : #include "VesTools.h" 25 : #include "VesBias.h" 26 : 27 : 28 : #include "tools/File.h" 29 : #include "core/ActionRegister.h" 30 : #include "core/PlumedMain.h" 31 : 32 : 33 : namespace PLMD { 34 : namespace ves { 35 : 36 : //+PLUMEDOC VES_UTILS VES_OUTPUT_FES 37 : /* 38 : Tool to output biases and free energy surfaces for VES biases from previously obtained coefficients. 39 : 40 : This action can be used to output to file biases and free energy surfaces for VES biases from 41 : previously obtained coefficients. It should be used through the \ref driver and 42 : can only be used in post processing. The VES bias needs to be defined in the 43 : exact same way as during the simulation. At the current moment this action does 44 : not support dynamic target distributions (e.g. well-tempered). 45 : 46 : \par Examples 47 : 48 : In the following input we define a VES bias and then read in the coefficient 49 : file coeffs.input.data and output the FES and bias every 500 iterations. 50 : 51 : \plumedfile 52 : phi: TORSION ATOMS=5,7,9,15 53 : psi: TORSION ATOMS=7,9,15,17 54 : 55 : bf1: BF_FOURIER ORDER=5 MINIMUM=-pi MAXIMUM=pi 56 : bf2: BF_FOURIER ORDER=5 MINIMUM=-pi MAXIMUM=pi 57 : 58 : VES_LINEAR_EXPANSION ... 59 : ARG=phi,psi 60 : BASIS_FUNCTIONS=bf1,bf2 61 : LABEL=ves1 62 : GRID_BINS=100,100 63 : PROJ_ARG1=phi 64 : PROJ_ARG2=psi 65 : ... VES_LINEAR_EXPANSION 66 : 67 : VES_OUTPUT_FES ... 68 : BIAS=ves1 69 : FES_OUTPUT=500 70 : FES_PROJ_OUTPUT=500 71 : BIAS_OUTPUT=500 72 : COEFFS_INPUT=coeffs.input.data 73 : ... VES_OUTPUT_FES 74 : \endplumedfile 75 : 76 : This input should be run through the driver by using a command similar to the 77 : following one where the trajectory/configuration file configuration.gro is needed to 78 : correctly define the CVs 79 : \verbatim 80 : plumed driver --plumed plumed.dat --igro configuration.gro 81 : \endverbatim 82 : 83 : */ 84 : //+ENDPLUMEDOC 85 : 86 3 : class OutputFesBias : public Action { 87 : 88 : public: 89 : static void registerKeywords(Keywords&); 90 : explicit OutputFesBias(const ActionOptions&); 91 0 : void update() {} 92 0 : void calculate() {} 93 0 : void apply() {} 94 : }; 95 : 96 : 97 7362 : PLUMED_REGISTER_ACTION(OutputFesBias,"VES_OUTPUT_FES") 98 : 99 : 100 4 : void OutputFesBias::registerKeywords(Keywords& keys) { 101 16 : keys.add("compulsory","BIAS","the label of the VES bias for to output the free energy surfaces and the bias files"); 102 16 : keys.add("compulsory","COEFFS_INPUT","the name of input coefficient file"); 103 16 : keys.add("optional","BIAS_OUTPUT","how often the bias(es) should be written out to file. Note that the value is given in terms of coefficient iterations."); 104 16 : keys.add("optional","FES_OUTPUT","how often the FES(s) should be written out to file. Note that the value is given in terms of coefficient iterations."); 105 16 : keys.add("optional","FES_PROJ_OUTPUT","how often the projections of the FES(s) should be written out to file. Note that the value is given in terms of coefficient iterations."); 106 : // 107 4 : } 108 : 109 : 110 3 : OutputFesBias::OutputFesBias(const ActionOptions&ao): 111 3 : Action(ao) 112 : { 113 : 114 3 : std::vector<std::string> bias_labels; 115 6 : parseVector("BIAS",bias_labels); 116 3 : if(bias_labels.size()>1) { 117 0 : plumed_merror(getName()+" only support one VES bias"); 118 : } 119 : 120 3 : std::string error_msg = ""; 121 6 : std::vector<VesBias*> bias_pntrs = VesTools::getPointersFromLabels<VesBias*>(bias_labels,plumed.getActionSet(),error_msg); 122 3 : if(error_msg.size()>0) {plumed_merror("Error in keyword BIAS of "+getName()+": "+error_msg);} 123 : 124 15 : for(unsigned int i=0; i<bias_pntrs.size(); i++) { 125 6 : if(bias_pntrs[i]->numberOfCoeffsSets()>1) { 126 0 : plumed_merror(getName()+" at the moment supports only VES biases with a single coefficient set"); 127 : } 128 : } 129 : 130 3 : std::vector<std::string> coeffs_fnames; 131 6 : parseVector("COEFFS_INPUT",coeffs_fnames); 132 3 : if(coeffs_fnames.size()!=bias_pntrs.size()) { 133 0 : plumed_merror(getName()+": there have to be as many coefficient file given in COEFFS_INPUT as VES biases given in BIAS"); 134 : } 135 : 136 3 : unsigned int bias_output_stride = 0; 137 6 : parse("BIAS_OUTPUT",bias_output_stride); 138 : 139 3 : unsigned int fes_output_stride = 0; 140 6 : parse("FES_OUTPUT",fes_output_stride); 141 : 142 3 : unsigned int fesproj_output_stride = 0; 143 6 : parse("FES_PROJ_OUTPUT",fesproj_output_stride); 144 : 145 3 : if(bias_output_stride == 0 && fes_output_stride == 0 && fesproj_output_stride == 0) { 146 0 : plumed_merror(getName()+": you are not telling the action to do anything, you need to use one of the keywords BIAS_OUTPUT, FES_OUTPUT, or FES_PROJ_OUTPUT"); 147 : } 148 : 149 3 : checkRead(); 150 : 151 15 : for(unsigned int i=0; i<bias_pntrs.size(); i++) { 152 : 153 6 : if(bias_pntrs[i]->dynamicTargetDistribution()) { 154 0 : plumed_merror(getName()+" does not support dynamic target distributions at the moment"); 155 : } 156 : 157 3 : if(bias_pntrs[i]->isStaticTargetDistFileOutputActive()) { 158 2 : bias_pntrs[i]->setupTargetDistFileOutput(); 159 2 : bias_pntrs[i]->writeTargetDistToFile(); 160 2 : bias_pntrs[i]->setupTargetDistProjFileOutput(); 161 2 : bias_pntrs[i]->writeTargetDistProjToFile(); 162 : } 163 : 164 : 165 3 : if(bias_output_stride>0) { 166 3 : bias_pntrs[i]->enableBiasFileOutput(); 167 3 : bias_pntrs[i]->setupBiasFileOutput(); 168 : } 169 : 170 3 : if(fes_output_stride>0) { 171 3 : bias_pntrs[i]->enableFesFileOutput(); 172 3 : bias_pntrs[i]->setupFesFileOutput(); 173 : } 174 : 175 3 : if(fesproj_output_stride>0) { 176 1 : bias_pntrs[i]->enableFesProjFileOutput(); 177 1 : bias_pntrs[i]->setupFesProjFileOutput(); 178 : } 179 : 180 3 : bias_pntrs[i]->enableIterationNumberInFilenames(); 181 : 182 6 : IFile ifile; 183 3 : ifile.open(coeffs_fnames[i]); 184 : 185 39 : while(ifile) { 186 : 187 36 : bias_pntrs[i]->resetBiasFileOutput(); 188 36 : bias_pntrs[i]->resetFesFileOutput(); 189 : 190 108 : if(bias_pntrs[i]->getCoeffsPntrs()[0]->readOneSetFromFile(ifile)>0) { 191 99 : unsigned int iteration = bias_pntrs[i]->getCoeffsPntrs()[0]->getIterationCounter(); 192 : 193 33 : if(bias_output_stride>0 && iteration%bias_output_stride==0) { 194 8 : bias_pntrs[i]->writeBiasToFile(); 195 : } 196 : 197 33 : if(fes_output_stride>0 && iteration%fes_output_stride==0) { 198 8 : bias_pntrs[i]->writeFesToFile(); 199 : } 200 : 201 33 : if(fesproj_output_stride>0 && iteration%fesproj_output_stride==0) { 202 3 : bias_pntrs[i]->writeFesProjToFile(); 203 : } 204 : 205 : } 206 : 207 : } 208 : 209 : } 210 : 211 3 : log.printf("Stopping"); 212 3 : plumed.stop(); 213 3 : } 214 : 215 : 216 : } 217 5517 : }