39int main(
int argc,
char **argv) {
41 return app.
Exec(argc, argv);
47 " options file for coarse graining")(
48 "trj-force", boost::program_options::value<string>(),
49 " coarse-grained trajectory containing forces of already known "
115 cout <<
"\nYou are using VOTCA!\n";
116 cout <<
"\nhey, somebody wants to forcematch!\n";
122 cout <<
"\nUsing constrained Least Squares!\n " << endl;
140 cout <<
"\nUsing simple Least Squares! " << endl;
163 throw runtime_error(
string(
"input format not supported: ") +
174 bool bonded_interaction,
181 bonded = bonded_interaction;
195 if (options->
exists(
"threebody")) {
204 if (options->
exists(
"fmatch.a")) {
205 a = options->
get(
"fmatch.a").
as<
double>();
207 if (options->
exists(
"fmatch.sigma")) {
208 sigma = options->
get(
"fmatch.sigma").
as<
double>();
210 if (options->
exists(
"fmatch.gamma")) {
211 gamma = options->
get(
"fmatch.gamma").
as<
double>();
222 if (options->
exists(
"fmatch.periodic")) {
228 std::cout <<
"a: " <<
a <<
" ,sigma: " <<
sigma <<
" ,gamma: " <<
gamma
232 double grid_min = options->
get(
"fmatch.min").
as<
double>();
233 double grid_max = options->
get(
"fmatch.max").
as<
double>();
234 double grid_step = options->
get(
"fmatch.step").
as<
double>();
241 cout <<
"Number of spline functions for the interaction " <<
splineName <<
":"
247 dx_out = options->
get(
"fmatch.out_step").
as<
double>();
270 cerr <<
"\nERROR in CGForceMatching::EndEvaluate - No blocks have been "
273 cerr <<
"It might be that you are using trajectory, which is smaller than "
274 "needed for one block"
276 cerr <<
"Check your input!" << endl;
280 cout <<
"\nWe are done, thank you very much!" << endl;
287 string file_extension =
".force";
288 string file_extension_pot =
".pot";
300 file_name = spline.splineName;
303 force_tab.
resize(spline.num_outgrid);
306 if (!(spline.threebody)) {
307 file_name = file_name + file_extension;
309 cout <<
"Updating file: " << file_name << endl;
315 if (spline.threebody) {
316 file_name = file_name + file_extension_pot;
317 file_nameDer = spline.splineName;
318 file_nameDer = file_nameDer + file_extension;
320 force_tabDer.
resize(spline.num_outgrid);
322 cout <<
"Updating files: " << file_name <<
" and: " << file_nameDer
326 spline.result = (spline.resSum).array() /
nblocks_;
327 spline.error = (((spline.resSum2).array() /
nblocks_ -
328 (spline.result).array().abs2()) /
333 if (spline.threebody) {
334 spline.resultDer = (spline.resSumDer).array() /
nblocks_;
335 spline.errorDer = (((spline.resSumDer2).array() /
nblocks_ -
336 (spline.resultDer).array().abs2()) /
343 double out_x = spline.Spline.getGridPoint(0);
348 if (!(spline.threebody)) {
350 force_tab.
set(i, out_x, (-1.0) * spline.result[i],
'i',
357 if (spline.threebody) {
359 force_tab.
set(i, out_x, (+1.0) * spline.result[i],
'i',
361 force_tabDer.
set(i, out_x, (-1.0) * spline.resultDer[i],
'i',
366 out_x += spline.dx_out;
369 force_tab.
Save(file_name);
375 if (spline.threebody) {
376 force_tabDer.
Save(file_nameDer);
378 force_tabDer.
clear();
385 throw std::runtime_error(
386 "CG Topology has 0 beads, check your mapping file!");
390 throw std::runtime_error(
391 "number of beads in topology and force topology does not match");
397 if (d.norm() >
dist_) {
399 throw std::runtime_error(
400 "One or more bead positions in mapped and reference force "
401 "trajectory differ by more than 1e-5");
410 if (sinfo.threebody) {
422 const Eigen::Vector3d &Force = conf->
getBead(iatom)->
getF();
430 throw std::runtime_error(
431 "\nERROR in csg_fmatch::EvalConfiguration - No forces in "
445 cout <<
"\nBlock No" <<
nblocks_ <<
" done!" << endl;
475 Eigen::HouseholderQR<Eigen::MatrixXd> dec(
A_);
477 Eigen::VectorXd residual =
b_ -
A_ *
x_;
480 double fm_resid = residual.cwiseAbs2().sum();
485 cout <<
"#### Force matching residual ####" << endl;
486 cout <<
" Chi_2[(kJ/(mol*nm))^2] = " << fm_resid << endl;
487 cout <<
"#################################" << endl;
498 sinfo.block_res_f =
x_.segment(mp, ngp);
499 sinfo.block_res_f2 =
x_.segment(mp + ngp, ngp);
502 sinfo.Spline.setSplineData(sinfo.block_res_f, sinfo.block_res_f2);
505 double out_x = sinfo.Spline.getGridPoint(0);
513 sinfo.resSum[i] += sinfo.Spline.Calculate(out_x);
516 sinfo.Spline.Calculate(out_x) * sinfo.Spline.Calculate(out_x);
520 if (sinfo.threebody) {
521 sinfo.resSumDer[i] += sinfo.Spline.CalculateDerivative(out_x);
523 sinfo.resSumDer2[i] += sinfo.Spline.CalculateDerivative(out_x) *
524 sinfo.Spline.CalculateDerivative(out_x);
528 if (i == grid_point_debug) {
529 cout <<
"This should be a number: " << sinfo.Spline.Calculate(out_x)
534 out_x += sinfo.dx_out;
550 sinfo.Spline.AddBCToFitMatrix(Matrix, line_tmp, col_tmp);
553 if (sinfo.periodic != 0) {
554 sinfo.Spline.AddBCSumZeroToFitMatrix(Matrix, line_tmp, col_tmp);
560 line_tmp += sfnum + 1;
561 col_tmp += 2 * (sfnum + 1);
573 std::vector<Interaction *> interList =
585 double var = inter->EvaluateVar(*conf);
588 for (
votca::Index loop = 0; loop < beads_in_int; loop++) {
590 Eigen::Vector3d gradient = inter->Grad(*conf, loop);
594 mpos, -gradient.x());
602 mpos, -gradient.z());
610 std::unique_ptr<NBList> nb;
612 bool gridsearch =
false;
617 }
else if (
options_.
get(
"cg.nbsearch").
as<
string>() ==
"simple") {
620 throw std::runtime_error(
"cg.nbsearch invalid, can be grid or simple");
624 nb = std::make_unique<NBListGrid>();
626 nb = std::make_unique<NBList>();
642 nb->Generate(beads1,
true);
644 nb->Generate(beads1, beads2,
true);
650 double var = pair->dist();
651 Eigen::Vector3d gradient = pair->r();
652 gradient.normalize();
674 mpos, -gradient.x());
682 mpos, -gradient.z());
691 std::unique_ptr<NBList_3Body> nb;
693 bool gridsearch =
false;
698 }
else if (
options_.
get(
"cg.nbsearch").
as<
string>() ==
"simple") {
701 throw std::runtime_error(
"cg.nbsearch invalid, can be grid or simple");
705 nb = std::make_unique<NBListGrid_3Body>();
707 nb = std::make_unique<NBList_3Body>();
710 nb->setCutoff(sinfo->
a);
727 nb->Generate(beads1,
true);
732 nb->Generate(beads1, beads2,
true);
740 nb->Generate(beads1, beads2, beads3,
true);
747 double distij = triple->dist12();
748 double distik = triple->dist13();
749 Eigen::Vector3d rij = triple->r12();
750 Eigen::Vector3d rik = triple->r13();
752 double gamma_sigma = (sinfo->
gamma) * (sinfo->
sigma);
753 double denomij = (distij - (sinfo->
a) * (sinfo->
sigma));
754 double denomik = (distik - (sinfo->
a) * (sinfo->
sigma));
755 double expij = std::exp(gamma_sigma / denomij);
756 double expik = std::exp(gamma_sigma / denomik);
763 std::acos(rij.dot(rik) / sqrt(rij.squaredNorm() * rik.squaredNorm()));
766 1.0 / (sqrt(1 - std::pow(rij.dot(rik), 2) /
767 (distij * distik * distij * distik)));
769 Eigen::Vector3d gradient1 =
771 ((rij + rik) / (distij * distik) -
772 rij.dot(rik) * (rik.squaredNorm() * rij + rij.squaredNorm() * rik) /
773 (distij * distij * distij * distik * distik * distik)) *
775 Eigen::Vector3d gradient2 =
776 ((rij / distij) * (gamma_sigma / (denomij * denomij)) +
777 (rik / distik) * (gamma_sigma / (denomik * denomik))) *
783 mpos, -gradient1.x(), -gradient2.x());
787 -gradient1.y(), -gradient2.y());
791 mpos, -gradient1.z(), -gradient2.z());
794 gradient1 = acos_prime *
795 (-rik / (distij * distik) +
796 rij.dot(rik) * rij / (distik * distij * distij * distij)) *
799 gradient2 = ((rij / distij) * (-1.0 * gamma_sigma / (denomij * denomij))) *
805 mpos, -gradient1.x(), -gradient2.x());
809 -gradient1.y(), -gradient2.y());
813 mpos, -gradient1.z(), -gradient2.z());
816 gradient1 = acos_prime *
817 (-rij / (distij * distik) +
818 rij.dot(rik) * rik / (distij * distik * distik * distik)) *
821 gradient2 = ((rik / distik) * (-1.0 * gamma_sigma / (denomik * denomik))) *
827 mpos, -gradient1.x(), -gradient2.x());
831 -gradient1.y(), -gradient2.y());
835 mpos, -gradient1.z(), -gradient2.z());
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
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.
void WriteOutFiles()
Write results to output files.
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.
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 Initialize(void) override
Initialize application data.
virtual const Eigen::Vector3d & getPos() const
Index Generate(Topology &top, const std::string &select)
Select all beads of type "select".
const Eigen::Vector3d & getF() const
get the force acting on the bead
bool HasF() const noexcept
bool EvaluateOptions() override
Process command line options.
void Initialize() override
Initialize application data.
base class for all interactions
topology of the whole system
std::vector< Interaction * > InteractionsInGroup(const std::string &group)
Bead * getBead(const Index i)
Returns a pointer to the bead with index i.
void CopyTopologyData(Topology *top)
copy topology data of different topology
int main(int argc, char **argv)
FileFormatFactory< TrajectoryReader > & TrjReaderFactory()
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
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)