votca 2024-dev
Loading...
Searching...
No Matches
Classes | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
CGForceMatching Class Reference

Implements force matching algorithm using cubic spline basis set. More...

#include <csg_fmatch.h>

Inheritance diagram for CGForceMatching:
Inheritance graph
[legend]
Collaboration diagram for CGForceMatching:
Collaboration graph
[legend]

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 DoMapping (void)
 overload and return true to enable mapping command line options
 
virtual bool DoMappingDefault (void)
 if DoMapping is true, will by default require mapping or not
 
virtual bool DoTrajectory (void)
 overload and return true to enable trajectory command line options
 
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)
 
virtual void BeginEvaluate (Topology *top, Topology *top_ref=nullptr)
 called before the first frame
 
virtual void EndEvaluate ()
 called after the last frame
 
virtual void EvalConfiguration (Topology *top, Topology *top_ref=nullptr)
 
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< WorkerForkWorker (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 ProgramName ()=0
 program name
 
virtual std::string VersionString ()
 version string of application
 
virtual void HelpText (std::ostream &out)=0
 help text of application without version information
 
virtual void Initialize ()=0
 Initialize application data.
 
virtual bool EvaluateOptions ()=0
 Process command line options.
 
void CheckRequired (const std::string &option_name, const std::string &error_msg="")
 Check weather required option is set.
 
virtual void Run ()=0
 Main body of application.
 
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 ()
 
virtual void ShowHelpText (std::ostream &out)
 

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< TrajectoryReadertrjreader_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< TrajectoryReadertraj_reader_
 
- Protected Attributes inherited from votca::tools::Application
std::map< std::string, boost::program_options::options_description > op_groups_
 
bool continue_execution_ = true
 

Detailed Description

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.

Todo:
force matching needs a big cleanup!

Definition at line 43 of file csg_fmatch.h.

Member Typedef Documentation

◆ SplineContainer

using CGForceMatching::SplineContainer = vector<SplineInfo>
protected

Definition at line 142 of file csg_fmatch.h.

Member Function Documentation

◆ BeginEvaluate()

void CGForceMatching::BeginEvaluate ( Topology top,
Topology top_atom 
)
overridevirtual

called before the first frame

Reimplemented from votca::csg::CsgApplication.

Definition at line 66 of file csg_fmatch.cc.

◆ DoMapping()

bool CGForceMatching::DoMapping ( void  )
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.

◆ DoTrajectory()

bool CGForceMatching::DoTrajectory ( void  )
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.

◆ EndEvaluate()

void CGForceMatching::EndEvaluate ( )
overridevirtual

called after the last frame

Reimplemented from votca::csg::CsgApplication.

Definition at line 267 of file csg_fmatch.cc.

◆ EvalBonded()

void CGForceMatching::EvalBonded ( Topology conf,
SplineInfo sinfo 
)
protected

For each trajectory frame writes equations for bonded interactions to matrix A_.

Definition at line 571 of file csg_fmatch.cc.

◆ EvalConfiguration()

void CGForceMatching::EvalConfiguration ( Topology conf,
Topology conf_atom = nullptr 
)
overridevirtual

called for each frame which is mapped

Reimplemented from votca::csg::CsgApplication.

Definition at line 383 of file csg_fmatch.cc.

◆ EvalNonbonded()

void CGForceMatching::EvalNonbonded ( Topology conf,
SplineInfo sinfo 
)
protected

For each trajectory frame writes equations for non-bonded interactions to matrix A_.

Definition at line 607 of file csg_fmatch.cc.

◆ EvalNonbonded_Threebody()

void CGForceMatching::EvalNonbonded_Threebody ( Topology conf,
SplineInfo sinfo 
)
protected

For each trajectory frame writes equations for non-bonded threebody interactions to matrix A_.

Definition at line 686 of file csg_fmatch.cc.

◆ EvaluateOptions()

bool CGForceMatching::EvaluateOptions ( )
overridevirtual

Process command line options.

Returns
true to continue, false to stop

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.

◆ FmatchAccumulateData()

void CGForceMatching::FmatchAccumulateData ( )
protected

Solves FM equations for one block and stores the results for further processing.

Definition at line 469 of file csg_fmatch.cc.

◆ FmatchAssignSmoothCondsToMatrix()

void CGForceMatching::FmatchAssignSmoothCondsToMatrix ( Eigen::MatrixXd &  Matrix)
protected

Assigns smoothing conditions to matrices A_ and B_constr_.

Definition at line 539 of file csg_fmatch.cc.

◆ HelpText()

void CGForceMatching::HelpText ( ostream &  out)
inlineoverridevirtual

help text of application without version information

Parameters
outostream for output

Implements votca::tools::Application.

Definition at line 46 of file csg_fmatch.h.

◆ Initialize()

void CGForceMatching::Initialize ( void  )
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.

◆ LoadOptions()

void CGForceMatching::LoadOptions ( const string &  file)

load options from the input file

Definition at line 565 of file csg_fmatch.cc.

◆ OpenForcesTrajectory()

void CGForceMatching::OpenForcesTrajectory ( )
protected

◆ ProgramName()

string CGForceMatching::ProgramName ( )
inlineoverridevirtual

program name

Returns
string with program name

overload this function to set the program name

Implements votca::tools::Application.

Definition at line 45 of file csg_fmatch.h.

◆ WriteOutFiles()

void CGForceMatching::WriteOutFiles ( )
protected

Write results to output files.

Definition at line 286 of file csg_fmatch.cc.

Member Data Documentation

◆ A_

Eigen::MatrixXd CGForceMatching::A_
protected

matrix used to store force matching equations

Definition at line 147 of file csg_fmatch.h.

◆ b_

Eigen::VectorXd CGForceMatching::b_
protected

vector used to store reference forces on CG beads (from atomistic simulations)

Definition at line 150 of file csg_fmatch.h.

◆ B_constr_

Eigen::MatrixXd CGForceMatching::B_constr_
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.

◆ bonded_

std::vector<votca::tools::Property *> CGForceMatching::bonded_
protected

list of bonded interactions

Definition at line 138 of file csg_fmatch.h.

◆ col_cntr_

votca::Index CGForceMatching::col_cntr_
protected

Definition at line 179 of file csg_fmatch.h.

◆ constr_least_sq_

bool CGForceMatching::constr_least_sq_
protected

Flag: true for constrained least squares, false for simple least squares.

Definition at line 169 of file csg_fmatch.h.

◆ dist_

double CGForceMatching::dist_
protected

accuracy for evaluating the difference in bead positions

Definition at line 165 of file csg_fmatch.h.

◆ frame_counter_

votca::Index CGForceMatching::frame_counter_
protected

Counter for trajectory frames.

Definition at line 160 of file csg_fmatch.h.

◆ has_existing_forces_

bool CGForceMatching::has_existing_forces_
protected

Definition at line 181 of file csg_fmatch.h.

◆ least_sq_offset_

votca::Index CGForceMatching::least_sq_offset_
protected

used in EvalConf to distinguish constrained and simple least squares

Definition at line 172 of file csg_fmatch.h.

◆ line_cntr_

votca::Index CGForceMatching::line_cntr_
protected

Counters for lines and columns in B_constr_.

Definition at line 179 of file csg_fmatch.h.

◆ nbeads_

votca::Index CGForceMatching::nbeads_
protected

Number of CG beads.

Definition at line 162 of file csg_fmatch.h.

◆ nblocks_

votca::Index CGForceMatching::nblocks_
protected

Current number of blocks.

Definition at line 176 of file csg_fmatch.h.

◆ nframes_

votca::Index CGForceMatching::nframes_
protected

Number of frames used in one block for block averaging.

Definition at line 174 of file csg_fmatch.h.

◆ nonbonded_

std::vector<votca::tools::Property *> CGForceMatching::nonbonded_
protected

list of non-bonded interactions

Definition at line 140 of file csg_fmatch.h.

◆ options_

votca::tools::Property CGForceMatching::options_
protected

Property object to handle input options.

Definition at line 136 of file csg_fmatch.h.

◆ splines_

SplineContainer CGForceMatching::splines_
protected

vector of SplineInfo * for all interactions

Definition at line 144 of file csg_fmatch.h.

◆ top_force_

Topology CGForceMatching::top_force_
protected

Definition at line 202 of file csg_fmatch.h.

◆ trjreader_force_

std::unique_ptr<TrajectoryReader> CGForceMatching::trjreader_force_
protected

Definition at line 203 of file csg_fmatch.h.

◆ x_

Eigen::VectorXd CGForceMatching::x_
protected

Solution of matrix equation A_ * x_ = b_ : CG force-field parameters.

Definition at line 153 of file csg_fmatch.h.


The documentation for this class was generated from the following files: