votca 2024.2-dev
|
Implements force matching algorithm using cubic spline basis set. More...
#include <csg_fmatch.h>
Classes | |
struct | SplineInfo |
structure, which contains CubicSpline object with related parameters More... | |
Public Member Functions | |
string | ProgramName () override |
program name | |
void | HelpText (ostream &out) override |
help text of application without version information | |
bool | DoTrajectory () override |
overload and return true to enable trajectory command line options | |
bool | DoMapping () override |
overload and return true to enable mapping command line options | |
void | Initialize (void) override |
Initialize application data. | |
bool | EvaluateOptions () override |
Process command line options. | |
void | BeginEvaluate (Topology *top, Topology *top_atom) override |
called before the first frame | |
void | EndEvaluate () override |
called after the last frame | |
void | EvalConfiguration (Topology *conf, Topology *conf_atom=nullptr) override |
called for each frame which is mapped | |
void | LoadOptions (const string &file) |
load options from the input file | |
Public Member Functions inherited from votca::csg::CsgApplication | |
CsgApplication ()=default | |
~CsgApplication () override=default | |
void | Initialize () override |
Initialize application data. | |
bool | EvaluateOptions () override |
Process command line options. | |
void | Run (void) override |
Main body of application. | |
void | ShowHelpText (std::ostream &out) override |
virtual bool | DoMappingDefault (void) |
if DoMapping is true, will by default require mapping or not | |
virtual bool | DoThreaded (void) |
virtual bool | SynchronizeThreads (void) |
virtual bool | NeedsTopology (void) |
if topology is always needed | |
virtual bool | EvaluateTopology (Topology *, Topology *=nullptr) |
called after topology was loaded | |
void | AddObserver (CGObserver *observer) |
bool | ProcessData (Worker *worker) |
Gets frames from TrajectoryReader in an ordered way and, if successful, calls Worker::EvalConfiguration for that frame. | |
virtual std::unique_ptr< Worker > | ForkWorker (void) |
virtual void | MergeWorker (Worker *worker) |
Public Member Functions inherited from votca::tools::Application | |
Application () | |
virtual | ~Application () |
int | Exec (int argc, char **argv) |
executes the program | |
virtual std::string | VersionString () |
version string of application | |
void | CheckRequired (const std::string &option_name, const std::string &error_msg="") |
Check weather required option is set. | |
boost::program_options::options_description_easy_init | AddProgramOptions (const std::string &group="") |
add option for command line | |
boost::program_options::variables_map & | OptionsMap () |
get available program options & descriptions | |
boost::program_options::options_description & | OptionsDesc () |
boost::program_options::options_description & | VisibleOptions () |
filters out the Hidden group from the options descriptions | |
void | StopExecution () |
call StopExecution after EvaluateOptions | |
Protected Types | |
using | SplineContainer = vector<SplineInfo> |
Protected Member Functions | |
void | FmatchAccumulateData () |
Solves FM equations for one block and stores the results for further processing. | |
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_. | |
void | EvalNonbonded (Topology *conf, SplineInfo *sinfo) |
For each trajectory frame writes equations for non-bonded interactions to matrix A_. | |
void | EvalNonbonded_Threebody (Topology *conf, SplineInfo *sinfo) |
For each trajectory frame writes equations for non-bonded threebody interactions to matrix A_. | |
void | WriteOutFiles () |
Write results to output files. | |
void | OpenForcesTrajectory () |
Protected Member Functions inherited from votca::tools::Application |
Protected Attributes | |
votca::tools::Property | options_ |
Property object to handle input options. | |
std::vector< votca::tools::Property * > | bonded_ |
list of bonded interactions | |
std::vector< votca::tools::Property * > | nonbonded_ |
list of non-bonded interactions | |
SplineContainer | splines_ |
vector of SplineInfo * for all interactions | |
Eigen::MatrixXd | A_ |
matrix used to store force matching equations | |
Eigen::VectorXd | b_ |
vector used to store reference forces on CG beads (from atomistic simulations) | |
Eigen::VectorXd | x_ |
Solution of matrix equation A_ * x_ = b_ : CG force-field parameters. | |
Eigen::MatrixXd | B_constr_ |
Additional matrix to handle constrained least squares fit contains constraints, which allow to get a real (smooth) spline (see VOTCA paper) | |
votca::Index | frame_counter_ |
Counter for trajectory frames. | |
votca::Index | nbeads_ |
Number of CG beads. | |
double | dist_ |
accuracy for evaluating the difference in bead positions | |
bool | constr_least_sq_ |
Flag: true for constrained least squares, false for simple least squares. | |
votca::Index | least_sq_offset_ |
used in EvalConf to distinguish constrained and simple least squares | |
votca::Index | nframes_ |
Number of frames used in one block for block averaging. | |
votca::Index | nblocks_ |
Current number of blocks. | |
votca::Index | line_cntr_ |
Counters for lines and columns in B_constr_. | |
votca::Index | col_cntr_ |
bool | has_existing_forces_ |
Topology | top_force_ |
std::unique_ptr< TrajectoryReader > | trjreader_force_ |
Protected Attributes inherited from votca::csg::CsgApplication | |
std::list< CGObserver * > | observers_ |
bool | do_mapping_ |
std::vector< std::unique_ptr< Worker > > | myWorkers_ |
Index | nframes_ |
bool | is_first_frame_ |
Index | nthreads_ |
tools::Mutex | nframesMutex_ |
tools::Mutex | traj_readerMutex_ |
std::vector< std::unique_ptr< tools::Mutex > > | threadsMutexesIn_ |
stores Mutexes used to impose order for input | |
std::vector< std::unique_ptr< tools::Mutex > > | threadsMutexesOut_ |
stores Mutexes used to impose order for output | |
std::unique_ptr< TrajectoryReader > | traj_reader_ |
Protected Attributes inherited from votca::tools::Application | |
std::map< std::string, boost::program_options::options_description > | op_groups_ |
bool | continue_execution_ = true |
Implements force matching algorithm using cubic spline basis set.
Force matching method to obtain a coarse-grained force field is implemented using cubic spline basis set. Block averaging over trajectory blocks is used for calculating CG forces and their errors.
Definition at line 43 of file csg_fmatch.h.
|
protected |
Definition at line 142 of file csg_fmatch.h.
called before the first frame
Reimplemented from votca::csg::CsgApplication.
Definition at line 66 of file csg_fmatch.cc.
|
inlineoverridevirtual |
overload and return true to enable mapping command line options
Reimplemented from votca::csg::CsgApplication.
Definition at line 51 of file csg_fmatch.h.
|
inlineoverridevirtual |
overload and return true to enable trajectory command line options
Reimplemented from votca::csg::CsgApplication.
Definition at line 50 of file csg_fmatch.h.
|
overridevirtual |
called after the last frame
Reimplemented from votca::csg::CsgApplication.
Definition at line 267 of file csg_fmatch.cc.
|
protected |
For each trajectory frame writes equations for bonded interactions to matrix A_.
Definition at line 571 of file csg_fmatch.cc.
|
overridevirtual |
called for each frame which is mapped
Reimplemented from votca::csg::CsgApplication.
Definition at line 383 of file csg_fmatch.cc.
|
protected |
For each trajectory frame writes equations for non-bonded interactions to matrix A_.
Definition at line 607 of file csg_fmatch.cc.
|
protected |
For each trajectory frame writes equations for non-bonded threebody interactions to matrix A_.
Definition at line 686 of file csg_fmatch.cc.
|
overridevirtual |
Process command line options.
EvaluateOptions is called by Run after parsing the command line. return true if everything is ok, false to stop and show help text.
Implements votca::tools::Application.
Definition at line 53 of file csg_fmatch.cc.
|
protected |
Solves FM equations for one block and stores the results for further processing.
Definition at line 469 of file csg_fmatch.cc.
|
protected |
Assigns smoothing conditions to matrices A_ and B_constr_.
Definition at line 539 of file csg_fmatch.cc.
|
inlineoverridevirtual |
help text of application without version information
out | ostream for output |
Implements votca::tools::Application.
Definition at line 46 of file csg_fmatch.h.
|
overridevirtual |
Initialize application data.
Initialize is called by run before parsing the command line. All necessary command line arguments can be added here
Implements votca::tools::Application.
Definition at line 44 of file csg_fmatch.cc.
void CGForceMatching::LoadOptions | ( | const string & | file | ) |
load options from the input file
Definition at line 565 of file csg_fmatch.cc.
|
protected |
|
inlineoverridevirtual |
program name
overload this function to set the program name
Implements votca::tools::Application.
Definition at line 45 of file csg_fmatch.h.
|
protected |
Write results to output files.
Definition at line 286 of file csg_fmatch.cc.
|
protected |
matrix used to store force matching equations
Definition at line 147 of file csg_fmatch.h.
|
protected |
vector used to store reference forces on CG beads (from atomistic simulations)
Definition at line 150 of file csg_fmatch.h.
|
protected |
Additional matrix to handle constrained least squares fit contains constraints, which allow to get a real (smooth) spline (see VOTCA paper)
Definition at line 157 of file csg_fmatch.h.
|
protected |
list of bonded interactions
Definition at line 138 of file csg_fmatch.h.
|
protected |
Definition at line 179 of file csg_fmatch.h.
|
protected |
Flag: true for constrained least squares, false for simple least squares.
Definition at line 169 of file csg_fmatch.h.
|
protected |
accuracy for evaluating the difference in bead positions
Definition at line 165 of file csg_fmatch.h.
|
protected |
Counter for trajectory frames.
Definition at line 160 of file csg_fmatch.h.
|
protected |
Definition at line 181 of file csg_fmatch.h.
|
protected |
used in EvalConf to distinguish constrained and simple least squares
Definition at line 172 of file csg_fmatch.h.
|
protected |
Counters for lines and columns in B_constr_.
Definition at line 179 of file csg_fmatch.h.
|
protected |
Number of CG beads.
Definition at line 162 of file csg_fmatch.h.
|
protected |
Current number of blocks.
Definition at line 176 of file csg_fmatch.h.
|
protected |
Number of frames used in one block for block averaging.
Definition at line 174 of file csg_fmatch.h.
|
protected |
list of non-bonded interactions
Definition at line 140 of file csg_fmatch.h.
|
protected |
Property object to handle input options.
Definition at line 136 of file csg_fmatch.h.
|
protected |
vector of SplineInfo * for all interactions
Definition at line 144 of file csg_fmatch.h.
|
protected |
Definition at line 202 of file csg_fmatch.h.
|
protected |
Definition at line 203 of file csg_fmatch.h.
|
protected |
Solution of matrix equation A_ * x_ = b_ : CG force-field parameters.
Definition at line 153 of file csg_fmatch.h.