18#ifndef VOTCA_CSG_CSG_FMATCH_H
19#define VOTCA_CSG_CSG_FMATCH_H
47 out <<
"Perform force matching (also called multiscale coarse-graining)";
62 Topology *conf_atom =
nullptr)
override;
138 std::vector<votca::tools::Property *>
bonded_;
Implements force matching algorithm using cubic spline basis set.
votca::Index nframes_
Number of frames used in one block for block averaging.
bool EvaluateOptions() override
Process command line options.
votca::Index nbeads_
Number of CG beads.
Eigen::VectorXd b_
vector used to store reference forces on CG beads (from atomistic simulations)
void FmatchAssignSmoothCondsToMatrix(Eigen::MatrixXd &Matrix)
Assigns smoothing conditions to matrices A_ and B_constr_.
void EvalBonded(Topology *conf, SplineInfo *sinfo)
For each trajectory frame writes equations for bonded interactions to matrix A_.
Eigen::MatrixXd A_
matrix used to store force matching equations
votca::Index frame_counter_
Counter for trajectory frames.
void EvalConfiguration(Topology *conf, Topology *conf_atom=nullptr) override
called for each frame which is mapped
votca::Index line_cntr_
Counters for lines and columns in B_constr_.
votca::tools::Property options_
Property object to handle input options.
Eigen::MatrixXd B_constr_
Additional matrix to handle constrained least squares fit contains constraints, which allow to get a ...
SplineContainer splines_
vector of SplineInfo * for all interactions
void BeginEvaluate(Topology *top, Topology *top_atom) override
called before the first frame
std::vector< votca::tools::Property * > nonbonded_
list of non-bonded interactions
string ProgramName() override
program name
bool has_existing_forces_
double dist_
accuracy for evaluating the difference in bead positions
votca::Index least_sq_offset_
used in EvalConf to distinguish constrained and simple least squares
void EvalNonbonded(Topology *conf, SplineInfo *sinfo)
For each trajectory frame writes equations for non-bonded interactions to matrix A_.
void FmatchAccumulateData()
Solves FM equations for one block and stores the results for further processing.
votca::Index nblocks_
Current number of blocks.
bool DoMapping() override
overload and return true to enable mapping command line options
void WriteOutFiles()
Write results to output files.
vector< SplineInfo > SplineContainer
void OpenForcesTrajectory()
std::vector< votca::tools::Property * > bonded_
list of bonded interactions
std::unique_ptr< TrajectoryReader > trjreader_force_
void EvalNonbonded_Threebody(Topology *conf, SplineInfo *sinfo)
For each trajectory frame writes equations for non-bonded threebody interactions to matrix A_.
bool constr_least_sq_
Flag: true for constrained least squares, false for simple least squares.
bool DoTrajectory() override
overload and return true to enable trajectory command line options
void EndEvaluate() override
called after the last frame
void LoadOptions(const string &file)
load options from the input file
Eigen::VectorXd x_
Solution of matrix equation A_ * x_ = b_ : CG force-field parameters.
void HelpText(ostream &out) override
help text of application without version information
void Initialize(void) override
Initialize application data.
topology of the whole system
structure, which contains CubicSpline object with related parameters
votca::Index num_outgrid
number of grid points for output
votca::Index num_gridpoints
number of spline grid points
votca::Index matr_pos
position in the A_ matrix (first coloumn which is occupied with this particular spline)
bool periodic
true if tabulated forces are periodic (e.g. for dihedral interactions)
votca::tools::Property * options_
pointer to Property object to handle input options
Eigen::VectorXd result
Final result: average over all blocks.
Eigen::VectorXd block_res_f
Result of 1 block calculation for f.
Eigen::VectorXd error
accuracy of the final result
string splineName
Spline Name.
double dx_out
dx for output. Calculated in the code
SplineInfo(votca::Index index, bool bonded_interaction, votca::Index matr_pos_col, votca::tools::Property *options)
constructor
Eigen::VectorXd resSum
sum of all block_res (used to calculate error)
Eigen::VectorXd resultDer
Final result of derivatives: average over all blocks.
votca::Index splineIndex
interaction index
pair< votca::Index, votca::Index > beadTypes
only for non-bonded interactions (seems like it is not used?)
bool threebody
true if non-bonded interaction is threebody interaction
Eigen::VectorXd errorDer
accuracy of the final result
Eigen::VectorXd block_res_f2
Result of 1 block calculation for f''.
bool bonded
true for bonded interactions, false for non-bonded
votca::tools::CubicSpline Spline
CubicSpline object.
double a
additional variables for treating cutoff of
Eigen::VectorXd resSumDer
sum of all block_res
votca::Index num_splinefun
number of spline functions
Eigen::VectorXd resSum2
sum of all squares of block_res (used to calculate error)
string type1
for non-bonded interactions: types of beads involved (type3 only used if threebody interaction)
Eigen::VectorXd resSumDer2
sum of all squares of block_res (used to calculate error)