LCOV - code coverage report
Current view: top level - ves - VesLinearExpansion.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 139 144 96.5 %
Date: 2024-10-18 14:00:25 Functions: 18 21 85.7 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             :    Copyright (c) 2016-2021 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 "VesBias.h"
      24             : #include "LinearBasisSetExpansion.h"
      25             : #include "CoeffsVector.h"
      26             : #include "CoeffsMatrix.h"
      27             : #include "BasisFunctions.h"
      28             : #include "Optimizer.h"
      29             : #include "TargetDistribution.h"
      30             : #include "VesTools.h"
      31             : 
      32             : #include "bias/Bias.h"
      33             : #include "core/ActionRegister.h"
      34             : #include "core/ActionSet.h"
      35             : #include "core/PlumedMain.h"
      36             : 
      37             : 
      38             : namespace PLMD {
      39             : namespace ves {
      40             : 
      41             : //+PLUMEDOC VES_BIAS VES_LINEAR_EXPANSION
      42             : /*
      43             : Linear basis set expansion bias.
      44             : 
      45             : This VES bias action takes the bias potential to be a linear expansion
      46             : in some basis set that is written as a product of one-dimensional basis functions.
      47             : For example, for one CV the bias would be written as
      48             : \f[
      49             : V(s_{1};\boldsymbol{\alpha}) = \sum_{i_{1}} \alpha_{i_{1}} \, f_{i_{1}}(s_{1}),
      50             : \f]
      51             : while for two CVs it is written as
      52             : \f[
      53             : V(s_{1},s_{2};\boldsymbol{\alpha}) = \sum_{i_{1},i_{2}} \alpha_{i_{1},i_{2}} \, f_{i_{1}}(s_{1}) \, f_{i_{2}}(s_{2})
      54             : \f]
      55             : where \f$\boldsymbol{\alpha}\f$ is the set of expansion coefficients that
      56             : are optimized within VES. With an appropriate choice of the basis functions
      57             : it is possible to represent any generic free energy surface.
      58             : The relationship between the bias and the free energy surface is given by
      59             : \f[
      60             : V(\mathbf{s}) = - F(\mathbf{s}) - \frac{1}{\beta} \log p(\mathbf{s}).
      61             : \f]
      62             : where \f$p(\mathbf{s})\f$ is the target distribution that is employed in the VES simulation.
      63             : 
      64             : \par Basis Functions
      65             : 
      66             : Various one-dimensional basis functions are available in the VES code,
      67             : see the complete list \ref ves_basisf "here".
      68             : At the current moment we recommend to use Legendre polynomials (\ref BF_LEGENDRE)
      69             : for non-periodic CVs and Fourier basis functions (\ref BF_FOURIER)
      70             : for periodic CV (e.g. dihedral angles).
      71             : 
      72             : To use basis functions within VES_LINEAR_EXPANSION you first need to
      73             : define them in the input file before the VES_LINEAR_EXPANSION action and
      74             : then give their labels using the BASIS_FUNCTIONS keyword.
      75             : 
      76             : 
      77             : \par Target Distributions
      78             : 
      79             : Various target distributions \f$p(\mathbf{s})\f$ are available in the VES code,
      80             : see the complete list \ref ves_targetdist "here".
      81             : 
      82             : To use a target distribution within VES_LINEAR_EXPANSION you first need to
      83             : define it in the input file before the VES_LINEAR_EXPANSION action and
      84             : then give its label using the TARGET_DISTRIBUTION keyword.
      85             : The default behavior if no TARGET_DISTRIBUTION is given is to
      86             : employ a uniform target distribution.
      87             : 
      88             : Some target distribution, like the well-tempered one (\ref TD_WELLTEMPERED),
      89             : are dynamic and need to be iteratively updated during the optimization.
      90             : 
      91             : \par Optimizer
      92             : 
      93             : In order to optimize the coefficients you will need to use VES_LINEAR_EXPANSION
      94             : in combination with an optimizer, see the list of optimizers available in the
      95             : VES code \ref ves_optimizer "here". At the current moment we recommend to
      96             : use the averaged stochastic gradient decent optimizer (\ref OPT_AVERAGED_SGD).
      97             : 
      98             : The optimizer should be defined after the VES_LINEAR_EXPANSION action.
      99             : 
     100             : \par Grid
     101             : 
     102             : Internally the code uses grids to calculate the basis set averages
     103             : over the target distribution that is needed for the gradient. The same grid is
     104             : also used for the output files (see next section).
     105             : The size of the grid is determined by the GRID_BINS keyword. By default it has
     106             : 100 grid points in each dimension, and generally this value should be sufficient.
     107             : 
     108             : \par Outputting Free Energy Surfaces and Other Files
     109             : 
     110             : It is possible to output on-the-fly during the simulation the free energy surface
     111             : estimated from the bias potential. How often this is done is specified within
     112             : the \ref ves_optimizer "optimizer" by using the FES_OUTPUT keyword. The filename
     113             : is specified by the FES_FILE keyword, but by default is it fes.LABEL.data,
     114             : with an added suffix indicating
     115             : the iteration number (iter-#).
     116             : 
     117             : For multi-dimensional case is it possible to also output projections of the
     118             : free energy surfaces. The arguments for which to do these projections is
     119             : specified using the numbered PROJ_ARG keywords. For these files a suffix
     120             : indicating the projection (proj-#) will be added to the filenames.
     121             : You will also need to specify the frequency of the output by using the
     122             : FES_PROJ_OUTPUT keyword within the optimizer.
     123             : 
     124             : It is also possible to output the bias potential itself, for this the relevant
     125             : keyword is BIAS_OUTPUT within the optimizer. The filename
     126             : is specified by the BIAS_FILE keyword, but by default is it bias.LABEL.data,
     127             : with an added suffix indicating the iteration number (iter-#).
     128             : 
     129             : Furthermore is it possible to output the target distribution, and its projections
     130             : (i.e. marginal distributions). The filenames of these files are specified with
     131             : the TARGETDIST_FILE, but by default is it targetdist.LABEL.data. The
     132             : logarithm of the target distribution will also be outputted to file that has the
     133             : added suffix log. For static target distribution these files will be outputted in
     134             : the beginning of the
     135             : simulation while for dynamic ones you will need to specify the frequency
     136             : of the output by using the TARGETDIST_OUTPUT and TARGETDIST_PROJ_OUTPUT
     137             : keywords within the optimizer.
     138             : 
     139             : It is also possible to output free energy surfaces and bias in post processing
     140             : by using the \ref VES_OUTPUT_FES action. However, be aware that this action
     141             : does does not support dynamic target distribution (e.g. well-tempered).
     142             : 
     143             : \par Static Bias
     144             : 
     145             : It is also possible to use VES_LINEAR_EXPANSION as a static bias that uses
     146             : previously obtained coefficients. In this case the coefficients should be
     147             : read in from the coefficient file given in the COEFFS keyword.
     148             : 
     149             : \par Bias Cutoff
     150             : 
     151             : It is possible to impose a cutoff on the bias potential using the procedure
     152             : introduced in \cite McCarty-PRL-2015 such that the free energy surface
     153             : is only flooded up to a certain value. The bias that results from this procedure
     154             : can then be used as a static bias for obtaining kinetic rates.
     155             : The value of the cutoff is given by the BIAS_CUTOFF keyword.
     156             : To impose the cutoff the code uses a Fermi switching function \f$1/(1+e^{\lambda x})\f$
     157             : where the parameter \f$\lambda\f$ controls how sharply the switchingfunction goes to zero.
     158             : The default value is \f$\lambda=10\f$ but this can be changed by using the
     159             : BIAS_CUTOFF_FERMI_LAMBDA keyword.
     160             : 
     161             : \par Examples
     162             : 
     163             : In the following example we run a VES_LINEAR_EXPANSION for one CV using
     164             : a Legendre basis functions (\ref BF_LEGENDRE) and a uniform target
     165             : distribution as no target distribution is specified. The coefficients
     166             : are optimized using averaged stochastic gradient descent optimizer
     167             : (\ref OPT_AVERAGED_SGD). Within the optimizer we specify that the
     168             : FES should be outputted to file every 500 coefficients iterations (the
     169             : FES_OUTPUT keyword).
     170             : Parameters that are very specific to the problem at hand, like the
     171             : order of the basis functions, the interval on which the
     172             : basis functions are defined, and the step size used
     173             : in the optimizer, are left unfilled.
     174             : \plumedfile
     175             : bf1: BF_LEGENDRE ORDER=__FILL__ MINIMUM=__FILL__ MAXIMUM=__FILL__
     176             : 
     177             : VES_LINEAR_EXPANSION ...
     178             :  ARG=d1
     179             :  BASIS_FUNCTIONS=bf1
     180             :  TEMP=__FILL__
     181             :  GRID_BINS=200
     182             :  LABEL=b1
     183             : ... VES_LINEAR_EXPANSION
     184             : 
     185             : OPT_AVERAGED_SGD ...
     186             :  BIAS=b1
     187             :  STRIDE=1000
     188             :  LABEL=o1
     189             :  STEPSIZE=__FILL__
     190             :  FES_OUTPUT=500
     191             :  COEFFS_OUTPUT=10
     192             : ... OPT_AVERAGED_SGD
     193             : \endplumedfile
     194             : 
     195             : In the following example we employ VES_LINEAR_EXPANSION for two CVs,
     196             : The first CV is periodic and therefore we employ a Fourier basis functions
     197             : (\ref BF_LEGENDRE) while the second CV is non-periodic so we employ a
     198             : Legendre polynomials as in the previous example. For the target distribution
     199             : we employ a well-tempered target distribution (\ref TD_WELLTEMPERED), which is
     200             : dynamic and needs to be iteratively updated with a stride that is given
     201             : using the TARGETDIST_STRIDE within the optimizer.
     202             : 
     203             : \plumedfile
     204             : bf1: BF_FOURIER  ORDER=__FILL__ MINIMUM=__FILL__ MAXIMUM=__FILL__
     205             : bf2: BF_LEGENDRE ORDER=__FILL__ MINIMUM=__FILL__ MAXIMUM=__FILL__
     206             : 
     207             : td_wt: TD_WELLTEMPERED BIASFACTOR=10.0
     208             : 
     209             : VES_LINEAR_EXPANSION ...
     210             :  ARG=cv1,cv2
     211             :  BASIS_FUNCTIONS=bf1,bf2
     212             :  TEMP=__FILL__
     213             :  GRID_BINS=100
     214             :  LABEL=b1
     215             :  TARGET_DISTRIBUTION=td_wt
     216             : ... VES_LINEAR_EXPANSION
     217             : 
     218             : OPT_AVERAGED_SGD ...
     219             :  BIAS=b1
     220             :  STRIDE=1000
     221             :  LABEL=o1
     222             :  STEPSIZE=__FILL__
     223             :  FES_OUTPUT=500
     224             :  COEFFS_OUTPUT=10
     225             :  TARGETDIST_STRIDE=500
     226             : ... OPT_AVERAGED_SGD
     227             : \endplumedfile
     228             : 
     229             : 
     230             : In the following example we employ a bias cutoff such that the bias
     231             : only fills the free energy landscape up a certain level. In this case
     232             : the target distribution is also dynamic and needs to iteratively updated.
     233             : 
     234             : \plumedfile
     235             : bf1: BF_LEGENDRE ORDER=__FILL__ MINIMUM=__FILL__ MAXIMUM=__FILL__
     236             : bf2: BF_LEGENDRE ORDER=__FILL__ MINIMUM=__FILL__ MAXIMUM=__FILL__
     237             : 
     238             : VES_LINEAR_EXPANSION ...
     239             :  ARG=cv1,cv2
     240             :  BASIS_FUNCTIONS=bf1,bf2
     241             :  TEMP=__FILL__
     242             :  GRID_BINS=100
     243             :  LABEL=b1
     244             :  BIAS_CUTOFF=20.0
     245             : ... VES_LINEAR_EXPANSION
     246             : 
     247             : OPT_AVERAGED_SGD ...
     248             :  BIAS=b1
     249             :  STRIDE=1000
     250             :  LABEL=o1
     251             :  STEPSIZE=__FILL__
     252             :  FES_OUTPUT=500
     253             :  COEFFS_OUTPUT=10
     254             :  TARGETDIST_STRIDE=500
     255             : ... OPT_AVERAGED_SGD
     256             : \endplumedfile
     257             : 
     258             : The optimized bias potential can then be used as a static bias for obtaining
     259             : kinetics. For this you need read in the final coefficients from file
     260             : (e.g. coeffs_final.data in this case) by using the
     261             : COEFFS keyword (also, no optimizer should be defined in the input)
     262             : \plumedfile
     263             : bf1: BF_LEGENDRE ORDER=__FILL__ MINIMUM=__FILL__ MAXIMUM=__FILL__
     264             : bf2: BF_LEGENDRE ORDER=__FILL__ MINIMUM=__FILL__ MAXIMUM=__FILL__
     265             : 
     266             : VES_LINEAR_EXPANSION ...
     267             :  ARG=cv1,cv2
     268             :  BASIS_FUNCTIONS=bf1,bf2
     269             :  TEMP=__FILL__
     270             :  GRID_BINS=100
     271             :  LABEL=b1
     272             :  BIAS_CUTOFF=20.0
     273             :  COEFFS=coeffs_final.data
     274             : ... VES_LINEAR_EXPANSION
     275             : \endplumedfile
     276             : 
     277             : 
     278             : 
     279             : */
     280             : //+ENDPLUMEDOC
     281             : 
     282             : 
     283             : class VesLinearExpansion : public VesBias {
     284             : private:
     285             :   unsigned int nargs_;
     286             :   std::vector<BasisFunctions*> basisf_pntrs_;
     287             :   std::unique_ptr<LinearBasisSetExpansion> bias_expansion_pntr_;
     288             :   size_t ncoeffs_;
     289             :   Value* valueForce2_;
     290             :   bool all_values_inside;
     291             :   std::vector<double> bf_values;
     292             :   bool bf_values_set;
     293             : public:
     294             :   explicit VesLinearExpansion(const ActionOptions&);
     295             :   ~VesLinearExpansion();
     296             :   void calculate() override;
     297             :   void update() override;
     298             :   void updateTargetDistributions() override;
     299             :   void restartTargetDistributions() override;
     300             :   //
     301             :   void setupBiasFileOutput() override;
     302             :   void writeBiasToFile() override;
     303             :   void resetBiasFileOutput() override;
     304             :   //
     305             :   void setupFesFileOutput() override;
     306             :   void writeFesToFile() override;
     307             :   void resetFesFileOutput() override;
     308             :   //
     309             :   void setupFesProjFileOutput() override;
     310             :   void writeFesProjToFile() override;
     311             :   //
     312             :   void writeTargetDistToFile() override;
     313             :   void writeTargetDistProjToFile() override;
     314             :   //
     315             :   double calculateReweightFactor() const override;
     316             :   //
     317             :   static void registerKeywords( Keywords& keys );
     318             : };
     319             : 
     320             : PLUMED_REGISTER_ACTION(VesLinearExpansion,"VES_LINEAR_EXPANSION")
     321             : 
     322          92 : void VesLinearExpansion::registerKeywords( Keywords& keys ) {
     323          92 :   VesBias::registerKeywords(keys);
     324             :   //
     325          92 :   VesBias::useInitialCoeffsKeywords(keys);
     326          92 :   VesBias::useTargetDistributionKeywords(keys);
     327          92 :   VesBias::useBiasCutoffKeywords(keys);
     328          92 :   VesBias::useGridBinKeywords(keys);
     329          92 :   VesBias::useProjectionArgKeywords(keys);
     330             :   //
     331          92 :   keys.use("ARG");
     332         184 :   keys.add("compulsory","BASIS_FUNCTIONS","the label of the one dimensional basis functions that should be used.");
     333         184 :   keys.addOutputComponent("force2","default","the instantaneous value of the squared force due to this bias potential.");
     334          92 : }
     335             : 
     336          90 : VesLinearExpansion::VesLinearExpansion(const ActionOptions&ao):
     337             :   PLUMED_VES_VESBIAS_INIT(ao),
     338          90 :   nargs_(getNumberOfArguments()),
     339          90 :   basisf_pntrs_(0),
     340          90 :   valueForce2_(NULL),
     341          90 :   all_values_inside(true),
     342          90 :   bf_values(0),
     343          90 :   bf_values_set(false)
     344             : {
     345             :   std::vector<std::string> basisf_labels;
     346          90 :   parseMultipleValues("BASIS_FUNCTIONS",basisf_labels,nargs_);
     347          90 :   checkRead();
     348             : 
     349          90 :   std::string error_msg = "";
     350         180 :   basisf_pntrs_ = VesTools::getPointersFromLabels<BasisFunctions*>(basisf_labels,plumed.getActionSet(),error_msg);
     351          90 :   if(error_msg.size()>0) {plumed_merror("Error in keyword BASIS_FUNCTIONS of "+getName()+": "+error_msg);}
     352             :   //
     353             : 
     354          90 :   std::vector<Value*> args_pntrs = getArguments();
     355             :   // check arguments and basis functions
     356             :   // this is done to avoid some issues with integration of target distribution
     357             :   // and periodic CVs, needs to be fixed later on.
     358         207 :   for(unsigned int i=0; i<args_pntrs.size(); i++) {
     359         117 :     if(args_pntrs[i]->isPeriodic() && !(basisf_pntrs_[i]->arePeriodic()) ) {
     360           0 :       plumed_merror("argument "+args_pntrs[i]->getName()+" is periodic while the basis functions " + basisf_pntrs_[i]->getLabel()+ " are not. You need to use the COMBINE action to remove the periodicity of the argument if you want to use these basis functions");
     361             :     }
     362         117 :     else if(!(args_pntrs[i]->isPeriodic()) && basisf_pntrs_[i]->arePeriodic() ) {
     363           1 :       log.printf("  warning: argument %s is not periodic while the basis functions %s used for it are periodic\n",args_pntrs[i]->getName().c_str(),basisf_pntrs_[i]->getLabel().c_str());
     364             :     }
     365             :   }
     366             : 
     367          90 :   addCoeffsSet(args_pntrs,basisf_pntrs_);
     368          90 :   ncoeffs_ = numberOfCoeffs();
     369          90 :   bool coeffs_read = readCoeffsFromFiles();
     370             : 
     371          90 :   checkThatTemperatureIsGiven();
     372         180 :   bias_expansion_pntr_ = Tools::make_unique<LinearBasisSetExpansion>(getLabel(),getBeta(),comm,args_pntrs,basisf_pntrs_,getCoeffsPntr());
     373          90 :   bias_expansion_pntr_->linkVesBias(this);
     374          90 :   bias_expansion_pntr_->setGridBins(this->getGridBins());
     375             :   //
     376          90 :   bf_values.assign(ncoeffs_,0.0);
     377             : 
     378             : 
     379             : 
     380          90 :   if(getNumberOfTargetDistributionPntrs()==0) {
     381          45 :     log.printf("  using an uniform target distribution: \n");
     382          45 :     bias_expansion_pntr_->setupUniformTargetDistribution();
     383             :     disableStaticTargetDistFileOutput();
     384             :   }
     385          45 :   else if(getNumberOfTargetDistributionPntrs()==1) {
     386          51 :     if(biasCutoffActive()) {getTargetDistributionPntrs()[0]->setupBiasCutoff();}
     387          45 :     bias_expansion_pntr_->setupTargetDistribution(getTargetDistributionPntrs()[0]);
     388          90 :     log.printf("  using target distribution of type %s with label %s \n",getTargetDistributionPntrs()[0]->getName().c_str(),getTargetDistributionPntrs()[0]->getLabel().c_str());
     389             :   }
     390             :   else {
     391           0 :     plumed_merror("problem with the TARGET_DISTRIBUTION keyword, either give no label or just one label.");
     392             :   }
     393          90 :   setTargetDistAverages(bias_expansion_pntr_->TargetDistAverages());
     394             :   //
     395          90 :   if(coeffs_read && biasCutoffActive()) {
     396           1 :     VesLinearExpansion::updateTargetDistributions();
     397             :   }
     398             :   //
     399          90 :   if(coeffs_read) {
     400           4 :     VesLinearExpansion::setupBiasFileOutput();
     401           4 :     VesLinearExpansion::writeBiasToFile();
     402             :   }
     403             : 
     404         180 :   addComponent("force2"); componentIsNotPeriodic("force2");
     405          90 :   valueForce2_=getPntrToComponent("force2");
     406          90 : }
     407             : 
     408             : 
     409         180 : VesLinearExpansion::~VesLinearExpansion() {
     410         270 : }
     411             : 
     412             : 
     413       23750 : void VesLinearExpansion::calculate() {
     414             : 
     415       23750 :   std::vector<double> cv_values(nargs_);
     416       23750 :   std::vector<double> forces(nargs_);
     417             : 
     418       60967 :   for(unsigned int k=0; k<nargs_; k++) {
     419       37217 :     cv_values[k]=getArgument(k);
     420             :   }
     421             : 
     422       23750 :   all_values_inside = true;
     423       23750 :   double bias = bias_expansion_pntr_->getBiasAndForces(cv_values,all_values_inside,forces,bf_values);
     424       23750 :   if(biasCutoffActive()) {
     425          63 :     applyBiasCutoff(bias,forces,bf_values);
     426          63 :     bf_values[0]=1.0;
     427             :   }
     428             :   double totalForce2 = 0.0;
     429       23750 :   if(all_values_inside) {
     430       60408 :     for(unsigned int k=0; k<nargs_; k++) {
     431       36852 :       setOutputForce(k,forces[k]);
     432       36852 :       totalForce2 += forces[k]*forces[k];
     433             :     }
     434             :   }
     435             : 
     436       23750 :   setBias(bias);
     437       23750 :   valueForce2_->set(totalForce2);
     438             : 
     439       23750 :   bf_values_set = true;
     440       23750 : }
     441             : 
     442             : 
     443       23750 : void VesLinearExpansion::update() {
     444       23750 :   if(!bf_values_set) {
     445           0 :     warning("VesLinearExpansion::update() is being called without calling VesLinearExpansion::calculate() first to calculate the basis function values. This can lead to incorrect behavior.");
     446             :   }
     447       23750 :   if(all_values_inside && bf_values_set) {
     448       23556 :     addToSampledAverages(bf_values);
     449             :   }
     450             :   std::fill(bf_values.begin(), bf_values.end(), 0.0);
     451       23750 :   bf_values_set = false;
     452       23750 : }
     453             : 
     454             : 
     455             : 
     456             : 
     457             : 
     458             : 
     459         355 : void VesLinearExpansion::updateTargetDistributions() {
     460         355 :   bias_expansion_pntr_->updateTargetDistribution();
     461         355 :   setTargetDistAverages(bias_expansion_pntr_->TargetDistAverages());
     462         355 : }
     463             : 
     464             : 
     465           8 : void VesLinearExpansion::restartTargetDistributions() {
     466          16 :   bias_expansion_pntr_->readInRestartTargetDistribution(getCurrentTargetDistOutputFilename());
     467           8 :   bias_expansion_pntr_->restartTargetDistribution();
     468           8 :   setTargetDistAverages(bias_expansion_pntr_->TargetDistAverages());
     469           8 : }
     470             : 
     471             : 
     472          87 : void VesLinearExpansion::setupBiasFileOutput() {
     473          87 :   bias_expansion_pntr_->setupBiasGrid(true);
     474          87 : }
     475             : 
     476             : 
     477         173 : void VesLinearExpansion::writeBiasToFile() {
     478         173 :   bias_expansion_pntr_->updateBiasGrid();
     479         519 :   auto ofile_pntr = getOFile(getCurrentBiasOutputFilename(),useMultipleWalkers());
     480         173 :   bias_expansion_pntr_->writeBiasGridToFile(*ofile_pntr);
     481         173 :   if(biasCutoffActive()) {
     482           5 :     bias_expansion_pntr_->updateBiasWithoutCutoffGrid();
     483          15 :     auto ofile_pntr2 = getOFile(getCurrentBiasOutputFilename("without-cutoff"),useMultipleWalkers());
     484           5 :     bias_expansion_pntr_->writeBiasWithoutCutoffGridToFile(*ofile_pntr2);
     485           5 :   }
     486         173 : }
     487             : 
     488             : 
     489          36 : void VesLinearExpansion::resetBiasFileOutput() {
     490             :   bias_expansion_pntr_->resetStepOfLastBiasGridUpdate();
     491          36 : }
     492             : 
     493             : 
     494          83 : void VesLinearExpansion::setupFesFileOutput() {
     495          83 :   bias_expansion_pntr_->setupFesGrid();
     496          83 : }
     497             : 
     498             : 
     499         169 : void VesLinearExpansion::writeFesToFile() {
     500         169 :   bias_expansion_pntr_->updateFesGrid();
     501         507 :   auto ofile_pntr = getOFile(getCurrentFesOutputFilename(),useMultipleWalkers());
     502         169 :   bias_expansion_pntr_->writeFesGridToFile(*ofile_pntr);
     503         169 : }
     504             : 
     505             : 
     506          36 : void VesLinearExpansion::resetFesFileOutput() {
     507             :   bias_expansion_pntr_->resetStepOfLastFesGridUpdate();
     508          36 : }
     509             : 
     510             : 
     511          17 : void VesLinearExpansion::setupFesProjFileOutput() {
     512          17 :   if(getNumberOfProjectionArguments()>0) {
     513           8 :     bias_expansion_pntr_->setupFesProjGrid();
     514             :   }
     515          17 : }
     516             : 
     517             : 
     518          36 : void VesLinearExpansion::writeFesProjToFile() {
     519          36 :   bias_expansion_pntr_->updateFesGrid();
     520          72 :   for(unsigned int i=0; i<getNumberOfProjectionArguments(); i++) {
     521             :     std::string suffix;
     522          36 :     Tools::convert(i+1,suffix);
     523          36 :     suffix = "proj-" + suffix;
     524          72 :     auto ofile_pntr = getOFile(getCurrentFesOutputFilename(suffix),useMultipleWalkers());
     525             :     std::vector<std::string> args = getProjectionArgument(i);
     526          36 :     bias_expansion_pntr_->writeFesProjGridToFile(args,*ofile_pntr);
     527          36 :   }
     528          36 : }
     529             : 
     530             : 
     531          82 : void VesLinearExpansion::writeTargetDistToFile() {
     532         164 :   auto ofile1_pntr = getOFile(getCurrentTargetDistOutputFilename(),useMultipleWalkers());
     533         164 :   auto ofile2_pntr = getOFile(getCurrentTargetDistOutputFilename("log"),useMultipleWalkers());
     534          82 :   bias_expansion_pntr_->writeTargetDistGridToFile(*ofile1_pntr);
     535          82 :   bias_expansion_pntr_->writeLogTargetDistGridToFile(*ofile2_pntr);
     536          82 : }
     537             : 
     538             : 
     539          13 : void VesLinearExpansion::writeTargetDistProjToFile() {
     540          33 :   for(unsigned int i=0; i<getNumberOfProjectionArguments(); i++) {
     541             :     std::string suffix;
     542          20 :     Tools::convert(i+1,suffix);
     543          20 :     suffix = "proj-" + suffix;
     544          40 :     auto ofile_pntr = getOFile(getCurrentTargetDistOutputFilename(suffix),useMultipleWalkers());
     545             :     std::vector<std::string> args = getProjectionArgument(i);
     546          20 :     bias_expansion_pntr_->writeTargetDistProjGridToFile(args,*ofile_pntr);
     547          20 :   }
     548          13 : }
     549             : 
     550             : 
     551           0 : double VesLinearExpansion::calculateReweightFactor() const {
     552           0 :   return bias_expansion_pntr_->calculateReweightFactor();
     553             : }
     554             : 
     555             : 
     556             : }
     557             : }

Generated by: LCOV version 1.16