votca 2024.2-dev
Loading...
Searching...
No Matches
csg_fmatch.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2021 The VOTCA Development Team (http://www.votca.org)
3 *
4 * Licensed under the Apache License, Version 2.0 (the "License");
5 * you may not use this file except in compliance with the License.
6 * You may obtain a copy of the License at
7 *
8 * http://www.apache.org/licenses/LICENSE-2.0
9 *
10 * Unless required by applicable law or agreed to in writing, software
11 * distributed under the License is distributed on an "AS IS" BASIS,
12 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13 * See the License for the specific language governing permissions and
14 * limitations under the License.
15 *
16 */
17
18#ifndef VOTCA_CSG_CSG_FMATCH_H
19#define VOTCA_CSG_CSG_FMATCH_H
20
21// VOTCA includes
24
25// Local VOTCA includes
28
29using namespace votca::csg;
30
31using namespace std;
32
44 public:
45 string ProgramName() override { return "csg_fmatch"; }
46 void HelpText(ostream &out) override {
47 out << "Perform force matching (also called multiscale coarse-graining)";
48 }
49
50 bool DoTrajectory() override { return true; }
51 bool DoMapping() override { return true; }
52
53 void Initialize(void) override;
54 bool EvaluateOptions() override;
55
57 void BeginEvaluate(Topology *top, Topology *top_atom) override;
59 void EndEvaluate() override;
61 void EvalConfiguration(Topology *conf,
62 Topology *conf_atom = nullptr) override;
64 void LoadOptions(const string &file);
65
66 protected:
69 struct SplineInfo {
71 SplineInfo(votca::Index index, bool bonded_interaction,
72 votca::Index matr_pos_col, votca::tools::Property *options);
82 bool bonded;
89 double a;
90 double sigma;
91 double gamma;
98 double dx_out;
100 pair<votca::Index, votca::Index> beadTypes;
101
103 Eigen::VectorXd block_res_f;
105 Eigen::VectorXd block_res_f2;
107 Eigen::VectorXd result;
109 Eigen::VectorXd error;
111 Eigen::VectorXd resSum;
113 Eigen::VectorXd resSum2;
114
115 // only needed for 3body nonbonded interactions as here force and potential
116 // are calculated simultaneously
118 Eigen::VectorXd resultDer;
120 Eigen::VectorXd errorDer;
122 Eigen::VectorXd resSumDer;
124 Eigen::VectorXd resSumDer2;
125
130 string type1, type2, type3; //
131
134 };
138 std::vector<votca::tools::Property *> bonded_;
140 std::vector<votca::tools::Property *> nonbonded_;
141
142 using SplineContainer = vector<SplineInfo>;
145
147 Eigen::MatrixXd A_;
150 Eigen::VectorXd b_;
153 Eigen::VectorXd x_; //
157 Eigen::MatrixXd B_constr_;
158
163
165 double dist_;
166
177
180
182
187 void FmatchAssignSmoothCondsToMatrix(Eigen::MatrixXd &Matrix);
190 void EvalBonded(Topology *conf, SplineInfo *sinfo);
193 void EvalNonbonded(Topology *conf, SplineInfo *sinfo);
196 void EvalNonbonded_Threebody(Topology *conf, SplineInfo *sinfo);
198 void WriteOutFiles();
199
201
203 std::unique_ptr<TrajectoryReader> trjreader_force_;
204};
205
206#endif // VOTCA_CSG_CSG_FMATCH_H
Implements force matching algorithm using cubic spline basis set.
Definition csg_fmatch.h:43
votca::Index nframes_
Number of frames used in one block for block averaging.
Definition csg_fmatch.h:174
bool EvaluateOptions() override
Process command line options.
Definition csg_fmatch.cc:53
votca::Index col_cntr_
Definition csg_fmatch.h:179
votca::Index nbeads_
Number of CG beads.
Definition csg_fmatch.h:162
Eigen::VectorXd b_
vector used to store reference forces on CG beads (from atomistic simulations)
Definition csg_fmatch.h:150
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
Definition csg_fmatch.h:147
votca::Index frame_counter_
Counter for trajectory frames.
Definition csg_fmatch.h:160
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_.
Definition csg_fmatch.h:179
votca::tools::Property options_
Property object to handle input options.
Definition csg_fmatch.h:136
Eigen::MatrixXd B_constr_
Additional matrix to handle constrained least squares fit contains constraints, which allow to get a ...
Definition csg_fmatch.h:157
SplineContainer splines_
vector of SplineInfo * for all interactions
Definition csg_fmatch.h:144
void BeginEvaluate(Topology *top, Topology *top_atom) override
called before the first frame
Definition csg_fmatch.cc:66
std::vector< votca::tools::Property * > nonbonded_
list of non-bonded interactions
Definition csg_fmatch.h:140
string ProgramName() override
program name
Definition csg_fmatch.h:45
bool has_existing_forces_
Definition csg_fmatch.h:181
double dist_
accuracy for evaluating the difference in bead positions
Definition csg_fmatch.h:165
votca::Index least_sq_offset_
used in EvalConf to distinguish constrained and simple least squares
Definition csg_fmatch.h:172
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.
Definition csg_fmatch.h:176
bool DoMapping() override
overload and return true to enable mapping command line options
Definition csg_fmatch.h:51
void WriteOutFiles()
Write results to output files.
vector< SplineInfo > SplineContainer
Definition csg_fmatch.h:142
void OpenForcesTrajectory()
Topology top_force_
Definition csg_fmatch.h:202
std::vector< votca::tools::Property * > bonded_
list of bonded interactions
Definition csg_fmatch.h:138
std::unique_ptr< TrajectoryReader > trjreader_force_
Definition csg_fmatch.h:203
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.
Definition csg_fmatch.h:169
bool DoTrajectory() override
overload and return true to enable trajectory command line options
Definition csg_fmatch.h:50
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.
Definition csg_fmatch.h:153
void HelpText(ostream &out) override
help text of application without version information
Definition csg_fmatch.h:46
void Initialize(void) override
Initialize application data.
Definition csg_fmatch.cc:44
topology of the whole system
Definition topology.h:81
A cubic spline class.
Definition cubicspline.h:54
class to manage program options with xml serialization functionality
Definition property.h:55
STL namespace.
Eigen::Index Index
Definition types.h:26
structure, which contains CubicSpline object with related parameters
Definition csg_fmatch.h:69
votca::Index num_outgrid
number of grid points for output
Definition csg_fmatch.h:78
votca::Index num_gridpoints
number of spline grid points
Definition csg_fmatch.h:76
votca::Index matr_pos
position in the A_ matrix (first coloumn which is occupied with this particular spline)
Definition csg_fmatch.h:96
bool periodic
true if tabulated forces are periodic (e.g. for dihedral interactions)
Definition csg_fmatch.h:85
votca::tools::Property * options_
pointer to Property object to handle input options
Definition csg_fmatch.h:133
Eigen::VectorXd result
Final result: average over all blocks.
Definition csg_fmatch.h:107
Eigen::VectorXd block_res_f
Result of 1 block calculation for f.
Definition csg_fmatch.h:103
Eigen::VectorXd error
accuracy of the final result
Definition csg_fmatch.h:109
string splineName
Spline Name.
Definition csg_fmatch.h:127
double dx_out
dx for output. Calculated in the code
Definition csg_fmatch.h:98
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)
Definition csg_fmatch.h:111
Eigen::VectorXd resultDer
Final result of derivatives: average over all blocks.
Definition csg_fmatch.h:118
votca::Index splineIndex
interaction index
Definition csg_fmatch.h:80
pair< votca::Index, votca::Index > beadTypes
only for non-bonded interactions (seems like it is not used?)
Definition csg_fmatch.h:100
bool threebody
true if non-bonded interaction is threebody interaction
Definition csg_fmatch.h:87
Eigen::VectorXd errorDer
accuracy of the final result
Definition csg_fmatch.h:120
Eigen::VectorXd block_res_f2
Result of 1 block calculation for f''.
Definition csg_fmatch.h:105
bool bonded
true for bonded interactions, false for non-bonded
Definition csg_fmatch.h:82
votca::tools::CubicSpline Spline
CubicSpline object.
Definition csg_fmatch.h:93
double a
additional variables for treating cutoff of
Definition csg_fmatch.h:89
Eigen::VectorXd resSumDer
sum of all block_res
Definition csg_fmatch.h:122
votca::Index num_splinefun
number of spline functions
Definition csg_fmatch.h:74
Eigen::VectorXd resSum2
sum of all squares of block_res (used to calculate error)
Definition csg_fmatch.h:113
string type1
for non-bonded interactions: types of beads involved (type3 only used if threebody interaction)
Definition csg_fmatch.h:130
Eigen::VectorXd resSumDer2
sum of all squares of block_res (used to calculate error)
Definition csg_fmatch.h:124