votca 2024-dev
Loading...
Searching...
No Matches
csg_fmatch.cc
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// Standard includes
19#include <cmath>
20#include <cstdio>
21#include <fstream>
22#include <iostream>
23#include <memory>
24#include <sstream>
25
26// VOTCA includes
28#include <votca/tools/linalg.h>
29#include <votca/tools/table.h>
30
31// Local VOTCA includes
32#include "votca/csg/beadlist.h"
35
36// Local private includes
37#include "csg_fmatch.h"
38
39int main(int argc, char **argv) {
41 return app.Exec(argc, argv);
42}
43
46 AddProgramOptions()("options", boost::program_options::value<string>(),
47 " options file for coarse graining")(
48 "trj-force", boost::program_options::value<string>(),
49 " coarse-grained trajectory containing forces of already known "
50 "interactions");
51}
52
55 CheckRequired("trj", "no trajectory file specified");
56 CheckRequired("options", "need to specify options file");
57 LoadOptions(OptionsMap()["options"].as<string>());
58
60 if (OptionsMap().count("trj-force")) {
62 }
63 return true;
64}
65
67 // set counters to zero value:
68 nblocks_ = 0;
70
71 // Number of CG beads in topology
72 nbeads_ = top->BeadCount();
73 // Set frame counter to zero
75
76 // accuracy for evaluating the difference in bead positions (default 1e-5)
77 dist_ = 1e-5;
78 if (options_.exists("cg.fmatch.dist")) {
79 dist_ = options_.get("cg.fmatch.dist").as<double>();
80 }
81
82 // read nframes_ from input file
83 nframes_ = options_.get("cg.fmatch.frames_per_block").as<votca::Index>();
84 // read constr_least_sq_ from input file
85 constr_least_sq_ = options_.get("cg.fmatch.constrainedLS").as<bool>();
86
87 // initializing bonded interactions
88 for (votca::tools::Property *prop : bonded_) {
89 // add spline to container
90 splines_.emplace_back(splines_.size(), true, col_cntr_, prop);
91 // adjust initial Eigen::Matrix3d dimensions:
92 line_cntr_ += splines_.back().num_gridpoints;
93 col_cntr_ += 2 * splines_.back().num_gridpoints;
94 // if periodic potential, one additional constraint has to be taken into
95 // account -> 1 additional line in matrix
96 if (splines_.back().periodic != 0) {
97 line_cntr_ += 1;
98 }
99 }
100
101 // initializing non-bonded interactions
102 for (votca::tools::Property *prop : nonbonded_) {
103 // add spline to container
104 splines_.emplace_back(splines_.size(), false, col_cntr_, prop);
105 // adjust initial Eigen::Matrix3d dimensions:
106 // number of constraints/restraints
107 line_cntr_ += splines_.back().num_gridpoints;
108 // number of coefficients
109 col_cntr_ += 2 * splines_.back().num_gridpoints;
110
111 // preliminary: use also spline functions for the threebody interaction. So
112 // far only angular interaction implemented
113 }
114
115 cout << "\nYou are using VOTCA!\n";
116 cout << "\nhey, somebody wants to forcematch!\n";
117
118 // now initialize A_, b_, x_ and probably B_constr_
119 // depending on least-squares algorithm used
120 if (constr_least_sq_) { // Constrained Least Squares
121
122 cout << "\nUsing constrained Least Squares!\n " << endl;
123
124 // assign least_sq_offset_
126
127 // resize and clear B_constr_
128 B_constr_ = Eigen::MatrixXd::Zero(line_cntr_, col_cntr_);
129
130 // resize Eigen::Matrix3d A_
131 A_ = Eigen::MatrixXd::Zero(3 * nbeads_ * nframes_, col_cntr_);
132 // resize vector b_
133 b_ = Eigen::VectorXd::Zero(3 * nbeads_ * nframes_);
134
135 // in case of constrained least squares smoothing conditions
136 // are assigned to Eigen::Matrix3d B_constr_
138 } else { // Simple Least Squares
139
140 cout << "\nUsing simple Least Squares! " << endl;
141 // assign least_sq_offset_
143
144 // resize Eigen::Matrix3d A_
145 A_ = Eigen::MatrixXd::Zero(line_cntr_ + 3 * nbeads_ * nframes_, col_cntr_);
146 // resize vector b_
147 b_ = Eigen::VectorXd::Zero(line_cntr_ + 3 * nbeads_ * nframes_);
148
149 // in case of simple least squares smoothing conditions
150 // are assigned to Eigen::Matrix3d A_
152 // clear b_ (only necessary in simple least squares)
153 b_.setZero();
154 }
155 // resize and clear x_
156 x_ = Eigen::VectorXd::Zero(col_cntr_);
157
161 TrjReaderFactory().Create(OptionsMap()["trj-force"].as<string>());
162 if (trjreader_force_ == nullptr) {
163 throw runtime_error(string("input format not supported: ") +
164 OptionsMap()["trj-force"].as<string>());
165 }
166 // open the trajectory
167 trjreader_force_->Open(OptionsMap()["trj-force"].as<string>());
168 // read in first frame
169 trjreader_force_->FirstFrame(top_force_);
170 }
171}
172
174 bool bonded_interaction,
175 votca::Index matr_pos_col,
176 votca::tools::Property *options) {
177 // initialize standard data
178 splineIndex = index;
179 options_ = options;
180 splineName = options->get("name").value();
181 bonded = bonded_interaction;
182 // in general natural boundary conditions are used for splines (default is no)
183 periodic = 0;
184 // check if non-bonded 3-body interaction or not (default is no)
185 threebody = 0;
186 // initialize additional parameters for threebody interactions
187 //(values of Molinero water potential)
188 a = 0.37; //(0.37 nm)
189 sigma = 1; //(dimensionless)
190 gamma = 0.12; //(0.12 nm = 1.2 Ang)
191
192 // get non-bonded information
193 if (!bonded) {
194 // check if option threebody exists
195 if (options->exists("threebody")) {
196 threebody = options->get("threebody").as<bool>();
197 }
198 // check if threebody interaction or not
199 if (threebody) {
200 type1 = options->get("type1").value();
201 type2 = options->get("type2").value();
202 type3 = options->get("type3").value();
203 // read in additional parameters for threebody interactions
204 if (options->exists("fmatch.a")) {
205 a = options->get("fmatch.a").as<double>();
206 }
207 if (options->exists("fmatch.sigma")) {
208 sigma = options->get("fmatch.sigma").as<double>();
209 }
210 if (options->exists("fmatch.gamma")) {
211 gamma = options->get("fmatch.gamma").as<double>();
212 }
213 }
214 // if not threebody only read in the two bead types
215 if (!threebody) {
216 type1 = options->get("type1").value();
217 type2 = options->get("type2").value();
218 }
219 }
220 if (bonded) {
221 // check if option periodic exists
222 if (options->exists("fmatch.periodic")) {
223 periodic = options->get("fmatch.periodic").as<bool>();
224 // set cubic spline Spline boundary conditions to periodic
225 Spline.setBCInt(1);
226 }
227 }
228 std::cout << "a: " << a << " ,sigma: " << sigma << " ,gamma: " << gamma
229 << std::endl;
230
231 // initialize the grid
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>();
235
236 // GenerateGrid returns number of grid points. We subtract 1 to get
237 // the number of spline functions
238 num_gridpoints = Spline.GenerateGrid(grid_min, grid_max, grid_step);
240
241 cout << "Number of spline functions for the interaction " << splineName << ":"
242 << num_splinefun << endl;
243
244 matr_pos = matr_pos_col;
245
246 // initialize grid for block averaging
247 dx_out = options->get("fmatch.out_step").as<double>();
248 // number of output grid points
249 num_outgrid = 1 + (votca::Index)((grid_max - grid_min) / dx_out);
250 result = Eigen::VectorXd::Zero(num_outgrid);
251 error = Eigen::VectorXd::Zero(num_outgrid);
252 resSum = Eigen::VectorXd::Zero(num_outgrid);
253 resSum2 = Eigen::VectorXd::Zero(num_outgrid);
254
255 // Only if threebody interaction, the derivatives are explicitly calculated
256 if (threebody) {
257 resultDer = Eigen::VectorXd::Zero(num_outgrid);
258 errorDer = Eigen::VectorXd::Zero(num_outgrid);
259 resSumDer = Eigen::VectorXd::Zero(num_outgrid);
260 resSumDer2 = Eigen::VectorXd::Zero(num_outgrid);
261 }
262
263 block_res_f = Eigen::VectorXd::Zero(num_outgrid);
264 block_res_f2 = Eigen::VectorXd::Zero(num_outgrid);
265}
266
268 // sanity check
269 if (nblocks_ == 0) {
270 cerr << "\nERROR in CGForceMatching::EndEvaluate - No blocks have been "
271 "processed so far"
272 << endl;
273 cerr << "It might be that you are using trajectory, which is smaller than "
274 "needed for one block"
275 << endl;
276 cerr << "Check your input!" << endl;
277 exit(-1);
278 }
279
280 cout << "\nWe are done, thank you very much!" << endl;
282 trjreader_force_->Close();
283 }
284}
285
287 string file_extension = ".force";
288 string file_extension_pot = ".pot";
289 string file_name;
290 string file_nameDer;
291 votca::tools::Table force_tab;
292 votca::tools::Table force_tabDer;
293
294 // table with error column
295 force_tab.SetHasYErr(true);
296 force_tabDer.SetHasYErr(true);
297
298 for (SplineInfo &spline : splines_) {
299 // construct meaningful outfile name
300 file_name = spline.splineName;
301
302 // resize table
303 force_tab.resize(spline.num_outgrid);
304
305 // If not threebody, the result represents the force
306 if (!(spline.threebody)) {
307 file_name = file_name + file_extension;
308 // print output file names on stdout
309 cout << "Updating file: " << file_name << endl;
310 }
311
312 // If threebody interaction, the result represents the potential and (-1)
313 // the derivative represents the force Only then, the derivatives are
314 // explicitly calculated
315 if (spline.threebody) {
316 file_name = file_name + file_extension_pot;
317 file_nameDer = spline.splineName;
318 file_nameDer = file_nameDer + file_extension;
319
320 force_tabDer.resize(spline.num_outgrid);
321 // print output file names on stdout
322 cout << "Updating files: " << file_name << " and: " << file_nameDer
323 << endl;
324 }
325
326 spline.result = (spline.resSum).array() / nblocks_;
327 spline.error = (((spline.resSum2).array() / nblocks_ -
328 (spline.result).array().abs2()) /
329 nblocks_)
330 .abs()
331 .sqrt();
332
333 if (spline.threebody) {
334 spline.resultDer = (spline.resSumDer).array() / nblocks_;
335 spline.errorDer = (((spline.resSumDer2).array() / nblocks_ -
336 (spline.resultDer).array().abs2()) /
337 nblocks_)
338 .abs()
339 .sqrt();
340 }
341
342 // first output point = first grid point
343 double out_x = spline.Spline.getGridPoint(0);
344 // loop over output grid
345 for (votca::Index i = 0; i < spline.num_outgrid; i++) {
346
347 // If not threebody the result is (-1) the force
348 if (!(spline.threebody)) {
349 // put point, result, flag and accuracy at point out_x into the table
350 force_tab.set(i, out_x, (-1.0) * spline.result[i], 'i',
351 spline.error[i]);
352 }
353
354 // If threebody interaction, force_tab represents the potential (-1) which
355 // is the Antiderivative of the force Only if threebody interaction, the
356 // derivatives are explicitly calculated
357 if (spline.threebody) {
358 // put point, result, flag and accuracy at point out_x into the table
359 force_tab.set(i, out_x, (+1.0) * spline.result[i], 'i',
360 spline.error[i]);
361 force_tabDer.set(i, out_x, (-1.0) * spline.resultDer[i], 'i',
362 spline.errorDer[i]);
363 }
364
365 // update out_x for the next iteration
366 out_x += spline.dx_out;
367 }
368 // save table in the file
369 force_tab.Save(file_name);
370
371 // clear the table for the next spline
372 force_tab.clear();
373
374 // Only if threebody interaction, the derivatives are explicitly calculated
375 if (spline.threebody) {
376 force_tabDer.Save(file_nameDer);
377 // clear the table for the next spline
378 force_tabDer.clear();
379 }
380 }
381}
382
384 if (conf->BeadCount() == 0) {
385 throw std::runtime_error(
386 "CG Topology has 0 beads, check your mapping file!");
387 }
389 if (conf->BeadCount() != top_force_.BeadCount()) {
390 throw std::runtime_error(
391 "number of beads in topology and force topology does not match");
392 }
393 for (votca::Index i = 0; i < conf->BeadCount(); ++i) {
394 conf->getBead(i)->F() -= top_force_.getBead(i)->getF();
395 Eigen::Vector3d d =
396 conf->getBead(i)->getPos() - top_force_.getBead(i)->getPos();
397 if (d.norm() > dist_) { // default is 1e-5, otherwise it can be a too
398 // strict criterion
399 throw std::runtime_error(
400 "One or more bead positions in mapped and reference force "
401 "trajectory differ by more than 1e-5");
402 }
403 }
404 }
405
406 for (SplineInfo &sinfo : splines_) {
407 if (sinfo.bonded) {
408 EvalBonded(conf, &sinfo);
409 } else {
410 if (sinfo.threebody) {
411 EvalNonbonded_Threebody(conf, &sinfo);
412 } else {
413 EvalNonbonded(conf, &sinfo);
414 }
415 }
416 }
417
418 // loop for the forces vector:
419 // hack, change the Has functions..
420 if (conf->getBead(0)->HasF()) {
421 for (votca::Index iatom = 0; iatom < nbeads_; ++iatom) {
422 const Eigen::Vector3d &Force = conf->getBead(iatom)->getF();
423 b_(least_sq_offset_ + 3 * nbeads_ * frame_counter_ + iatom) = Force.x();
425 Force.y();
427 iatom) = Force.z();
428 }
429 } else {
430 throw std::runtime_error(
431 "\nERROR in csg_fmatch::EvalConfiguration - No forces in "
432 "configuration!");
433 }
434 // update the frame counter
435 frame_counter_ += 1;
436
437 if (frame_counter_ % nframes_ == 0) { // at this point we processed nframes_
438 // frames, which is enough for one
439 // block
440 // update block counter
441 nblocks_++;
442 // solve FM equations and accumulate the result
444 // print status information
445 cout << "\nBlock No" << nblocks_ << " done!" << endl;
446 // write results to output files
448
449 // we must count frames from zero again for the next block
450 frame_counter_ = 0;
451 if (constr_least_sq_) { // Constrained Least Squares
452 // Matrices should be cleaned after each block is evaluated
453 A_.setZero();
454 b_.setZero();
455 // clear and assign smoothing conditions to B_constr_
457 } else { // Simple Least Squares
458 // Matrices should be cleaned after each block is evaluated
459 // clear and assign smoothing conditions to A_
461 b_.setZero();
462 }
463 }
465 trjreader_force_->NextFrame(top_force_);
466 }
467}
468
470 if (constr_least_sq_) { // Constrained Least Squares
471 // Solving linear equations system
473 } else { // Simple Least Squares
474
475 Eigen::HouseholderQR<Eigen::MatrixXd> dec(A_);
476 x_ = dec.solve(b_);
477 Eigen::VectorXd residual = b_ - A_ * x_;
478 // calculate FM residual - quality of FM
479 // FM residual is calculated in (kJ/(mol*nm))^2
480 double fm_resid = residual.cwiseAbs2().sum();
481
482 fm_resid /= (double)(3 * nbeads_ * frame_counter_);
483
484 cout << endl;
485 cout << "#### Force matching residual ####" << endl;
486 cout << " Chi_2[(kJ/(mol*nm))^2] = " << fm_resid << endl;
487 cout << "#################################" << endl;
488 cout << endl;
489 }
490
491 for (SplineInfo &sinfo : splines_) {
492 votca::Index mp = sinfo.matr_pos;
493 votca::Index ngp = sinfo.num_gridpoints;
494
495 // x_ contains results for all splines. Here we cut the results for one
496 // spline
497
498 sinfo.block_res_f = x_.segment(mp, ngp);
499 sinfo.block_res_f2 = x_.segment(mp + ngp, ngp);
500
501 // result cut before is assigned to the corresponding spline
502 sinfo.Spline.setSplineData(sinfo.block_res_f, sinfo.block_res_f2);
503
504 // first output point = first grid point
505 double out_x = sinfo.Spline.getGridPoint(0);
506
507 // point in the middle of the output grid for printing debug information
508 votca::Index grid_point_debug = sinfo.num_outgrid / 2;
509
510 // loop over output grid
511 for (votca::Index i = 0; i < sinfo.num_outgrid; i++) {
512 // update resSum (add result of a particular block)
513 sinfo.resSum[i] += sinfo.Spline.Calculate(out_x);
514 // update resSum2 (add result of a particular block)
515 sinfo.resSum2[i] +=
516 sinfo.Spline.Calculate(out_x) * sinfo.Spline.Calculate(out_x);
517
518 // Only if threebody interaction, the derivatives are explicitly
519 // calculated
520 if (sinfo.threebody) {
521 sinfo.resSumDer[i] += sinfo.Spline.CalculateDerivative(out_x);
522 // update resSumDer2 (add result of a particular block)
523 sinfo.resSumDer2[i] += sinfo.Spline.CalculateDerivative(out_x) *
524 sinfo.Spline.CalculateDerivative(out_x);
525 }
526
527 // print useful debug information
528 if (i == grid_point_debug) {
529 cout << "This should be a number: " << sinfo.Spline.Calculate(out_x)
530 << " " << endl;
531 }
532
533 // output point for the next iteration
534 out_x += sinfo.dx_out;
535 }
536 }
537}
538
540 // This function assigns Spline smoothing conditions to the Matrix.
541 // For the simple least squares the function is used for Eigen::Matrix3d A_
542 // For constrained least squares - for Eigen::Matrix3d B_constr_
543
544 Matrix.setZero();
545 votca::Index line_tmp = 0;
546 votca::Index col_tmp = 0;
547
548 for (SplineInfo &sinfo : splines_) {
549
550 sinfo.Spline.AddBCToFitMatrix(Matrix, line_tmp, col_tmp);
551 // if periodic potential, one additional constraint has to be taken into
552 // account!
553 if (sinfo.periodic != 0) {
554 sinfo.Spline.AddBCSumZeroToFitMatrix(Matrix, line_tmp, col_tmp);
555 // update counter
556 line_tmp += 1;
557 }
558 // update counters
559 votca::Index sfnum = sinfo.num_splinefun;
560 line_tmp += sfnum + 1;
561 col_tmp += 2 * (sfnum + 1);
562 }
563}
564
565void CGForceMatching::LoadOptions(const string &file) {
566 options_.LoadFromXML(file);
567 bonded_ = options_.Select("cg.bonded");
568 nonbonded_ = options_.Select("cg.non-bonded");
569}
570
572
573 std::vector<Interaction *> interList =
574 conf->InteractionsInGroup(sinfo->splineName);
575
576 for (Interaction *inter : interList) {
577
578 votca::Index beads_in_int = inter->BeadCount(); // 2 for bonds, 3 for
579 // angles, 4 for dihedrals
580
582
583 votca::Index mpos = sinfo->matr_pos;
584
585 double var = inter->EvaluateVar(*conf); // value of bond, angle,
586 // or dihedral
587
588 for (votca::Index loop = 0; loop < beads_in_int; loop++) {
589 votca::Index ii = inter->getBeadId(loop);
590 Eigen::Vector3d gradient = inter->Grad(*conf, loop);
591
592 SP.AddToFitMatrix(A_, var,
594 mpos, -gradient.x());
596 A_, var,
597 least_sq_offset_ + 3 * nbeads_ * frame_counter_ + nbeads_ + ii, mpos,
598 -gradient.y());
600 A_, var,
602 mpos, -gradient.z());
603 }
604 }
605}
606
608
609 // generate the neighbour list
610 std::unique_ptr<NBList> nb;
611
612 bool gridsearch = false;
613
614 if (options_.exists("cg.nbsearch")) {
615 if (options_.get("cg.nbsearch").as<string>() == "grid") {
616 gridsearch = true;
617 } else if (options_.get("cg.nbsearch").as<string>() == "simple") {
618 gridsearch = false;
619 } else {
620 throw std::runtime_error("cg.nbsearch invalid, can be grid or simple");
621 }
622 }
623 if (gridsearch) {
624 nb = std::make_unique<NBListGrid>();
625 } else {
626 nb = std::make_unique<NBList>();
627 }
628
629 nb->setCutoff(
630 sinfo->options_->get("fmatch.max").as<double>()); // implement different
631 // cutoffs for
632 // different
633 // interactions!
634
635 // generate the bead lists
636 BeadList beads1, beads2;
637 beads1.Generate(*conf, sinfo->type1);
638 beads2.Generate(*conf, sinfo->type2);
639
640 // is it same types or different types?
641 if (sinfo->type1 == sinfo->type2) {
642 nb->Generate(beads1, true);
643 } else {
644 nb->Generate(beads1, beads2, true);
645 }
646
647 for (BeadPair *pair : *nb) {
648 votca::Index iatom = pair->first()->getId();
649 votca::Index jatom = pair->second()->getId();
650 double var = pair->dist();
651 Eigen::Vector3d gradient = pair->r();
652 gradient.normalize();
653
655
656 votca::Index mpos = sinfo->matr_pos;
657
658 // add iatom
659 SP.AddToFitMatrix(A_, var,
661 mpos, gradient.x());
663 A_, var,
664 least_sq_offset_ + 3 * nbeads_ * frame_counter_ + nbeads_ + iatom, mpos,
665 gradient.y());
667 A_, var,
668 least_sq_offset_ + 3 * nbeads_ * frame_counter_ + 2 * nbeads_ + iatom,
669 mpos, gradient.z());
670
671 // add jatom
672 SP.AddToFitMatrix(A_, var,
674 mpos, -gradient.x());
676 A_, var,
677 least_sq_offset_ + 3 * nbeads_ * frame_counter_ + nbeads_ + jatom, mpos,
678 -gradient.y());
680 A_, var,
681 least_sq_offset_ + 3 * nbeads_ * frame_counter_ + 2 * nbeads_ + jatom,
682 mpos, -gradient.z());
683 }
684}
685
687 SplineInfo *sinfo) {
688 // so far option gridsearch ignored. Only simple search
689
690 // generate the neighbour list
691 std::unique_ptr<NBList_3Body> nb;
692
693 bool gridsearch = false;
694
695 if (options_.exists("cg.nbsearch")) {
696 if (options_.get("cg.nbsearch").as<string>() == "grid") {
697 gridsearch = true;
698 } else if (options_.get("cg.nbsearch").as<string>() == "simple") {
699 gridsearch = false;
700 } else {
701 throw std::runtime_error("cg.nbsearch invalid, can be grid or simple");
702 }
703 }
704 if (gridsearch) {
705 nb = std::make_unique<NBListGrid_3Body>();
706 } else {
707 nb = std::make_unique<NBList_3Body>();
708 }
709
710 nb->setCutoff(sinfo->a); // implement different cutoffs for different
711 // interactions!
712 // Here, a is the distance between two beads of a triple, where the 3-body
713 // interaction is zero
714
715 // generate the bead lists
716 BeadList beads1, beads2, beads3;
717 beads1.Generate(*conf, sinfo->type1);
718 beads2.Generate(*conf, sinfo->type2);
719 beads3.Generate(*conf, sinfo->type3);
720
721 // Generate the 3body neighbour lists
722 // check if type2 and type3 are the same
723 if (sinfo->type2 == sinfo->type3) {
724 // if then type2 and type1 are the same, all three types are the same
725 // use the Generate function for this case
726 if (sinfo->type1 == sinfo->type2) {
727 nb->Generate(beads1, true);
728 }
729 // else use the Generate function for type2 being equal to type3 (and type1
730 // being different)
731 if (sinfo->type1 != sinfo->type2) {
732 nb->Generate(beads1, beads2, true);
733 }
734 }
735 // If type2 and type3 are not the same, use the Generate function for three
736 // different bead types (Even if type1 and type2 or type1 and type3 are the
737 // same, the Generate function for two different beadtypes is only applicable
738 // for the case that type2 is equal to type3
739 if (sinfo->type2 != sinfo->type3) {
740 nb->Generate(beads1, beads2, beads3, true);
741 }
742
743 for (BeadTriple *triple : *nb) {
744 votca::Index iatom = triple->bead1()->getId();
745 votca::Index jatom = triple->bead2()->getId();
746 votca::Index katom = triple->bead3()->getId();
747 double distij = triple->dist12();
748 double distik = triple->dist13();
749 Eigen::Vector3d rij = triple->r12();
750 Eigen::Vector3d rik = triple->r13();
751
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);
757
759
760 votca::Index mpos = sinfo->matr_pos;
761
762 double var =
763 std::acos(rij.dot(rik) / sqrt(rij.squaredNorm() * rik.squaredNorm()));
764
765 double acos_prime =
766 1.0 / (sqrt(1 - std::pow(rij.dot(rik), 2) /
767 (distij * distik * distij * distik)));
768
769 Eigen::Vector3d gradient1 =
770 acos_prime *
771 ((rij + rik) / (distij * distik) -
772 rij.dot(rik) * (rik.squaredNorm() * rij + rij.squaredNorm() * rik) /
773 (distij * distij * distij * distik * distik * distik)) *
774 expij * expik;
775 Eigen::Vector3d gradient2 =
776 ((rij / distij) * (gamma_sigma / (denomij * denomij)) +
777 (rik / distik) * (gamma_sigma / (denomik * denomik))) *
778 expij * expik;
779
780 // add iatom
781 SP.AddToFitMatrix(A_, var,
783 mpos, -gradient1.x(), -gradient2.x());
785 A_, var,
786 least_sq_offset_ + 3 * nbeads_ * frame_counter_ + nbeads_ + iatom, mpos,
787 -gradient1.y(), -gradient2.y());
789 A_, var,
790 least_sq_offset_ + 3 * nbeads_ * frame_counter_ + 2 * nbeads_ + iatom,
791 mpos, -gradient1.z(), -gradient2.z());
792
793 // evaluate gradient1 and gradient2 for jatom:
794 gradient1 = acos_prime *
795 (-rik / (distij * distik) +
796 rij.dot(rik) * rij / (distik * distij * distij * distij)) *
797 expij * expik;
798 // gradient2
799 gradient2 = ((rij / distij) * (-1.0 * gamma_sigma / (denomij * denomij))) *
800 expij * expik;
801
802 // add jatom
803 SP.AddToFitMatrix(A_, var,
805 mpos, -gradient1.x(), -gradient2.x());
807 A_, var,
808 least_sq_offset_ + 3 * nbeads_ * frame_counter_ + nbeads_ + jatom, mpos,
809 -gradient1.y(), -gradient2.y());
811 A_, var,
812 least_sq_offset_ + 3 * nbeads_ * frame_counter_ + 2 * nbeads_ + jatom,
813 mpos, -gradient1.z(), -gradient2.z());
814
815 // evaluate gradient1 and gradient2 for katom:
816 gradient1 = acos_prime *
817 (-rij / (distij * distik) +
818 rij.dot(rik) * rik / (distij * distik * distik * distik)) *
819 expij * expik;
820 // gradient2
821 gradient2 = ((rik / distik) * (-1.0 * gamma_sigma / (denomik * denomik))) *
822 expij * expik;
823
824 // add katom
825 SP.AddToFitMatrix(A_, var,
827 mpos, -gradient1.x(), -gradient2.x());
829 A_, var,
830 least_sq_offset_ + 3 * nbeads_ * frame_counter_ + nbeads_ + katom, mpos,
831 -gradient1.y(), -gradient2.y());
833 A_, var,
834 least_sq_offset_ + 3 * nbeads_ * frame_counter_ + 2 * nbeads_ + katom,
835 mpos, -gradient1.z(), -gradient2.z());
836 }
837}
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
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
void WriteOutFiles()
Write results to output files.
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
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 Initialize(void) override
Initialize application data.
Definition csg_fmatch.cc:44
virtual const Eigen::Vector3d & getPos() const
Definition basebead.h:166
Index Generate(Topology &top, const std::string &select)
Select all beads of type "select".
Definition beadlist.cc:30
A particle pair.
Definition beadpair.h:35
A triplet of tree Beads.
Definition beadtriple.h:32
const Eigen::Vector3d & getF() const
get the force acting on the bead
Definition bead.h:361
bool HasF() const noexcept
Definition bead.h:235
Eigen::Vector3d & F()
Definition bead.h:210
bool EvaluateOptions() override
Process command line options.
void Initialize() override
Initialize application data.
base class for all interactions
Definition interaction.h:40
topology of the whole system
Definition topology.h:81
std::vector< Interaction * > InteractionsInGroup(const std::string &group)
Definition topology.cc:201
Index BeadCount() const
Definition topology.h:150
Bead * getBead(const Index i)
Returns a pointer to the bead with index i.
Definition topology.h:227
void CopyTopologyData(Topology *top)
copy topology data of different topology
Definition topology.cc:102
int Exec(int argc, char **argv)
executes the program
boost::program_options::variables_map & OptionsMap()
get available program options & descriptions
boost::program_options::options_description_easy_init AddProgramOptions(const std::string &group="")
add option for command line
void CheckRequired(const std::string &option_name, const std::string &error_msg="")
Check weather required option is set.
A cubic spline class.
Definition cubicspline.h:54
void AddToFitMatrix(matrix_type &M, double x, Index offset1, Index offset2=0, double scale=1)
Add a point (one entry) to fitting matrix.
class to manage program options with xml serialization functionality
Definition property.h:55
Property & get(const std::string &key)
get existing property
Definition property.cc:79
std::string & value()
reference to value of property
Definition property.h:153
bool exists(const std::string &key) const
check whether property exists
Definition property.cc:122
T as() const
return value as type
Definition property.h:283
std::vector< Property * > Select(const std::string &filter)
select property based on a filter
Definition property.cc:185
void LoadFromXML(std::string filename)
Definition property.cc:238
Spline Class.
Definition spline.h:36
void setBCInt(Index bc)
Set the boundary type of the spline.
Definition spline.h:89
Index GenerateGrid(double min, double max, double h)
Generate the grid for fitting from "min" to "max" in steps of "h".
Definition spline.cc:26
class to store tables like rdfs, tabulated potentials, etc
Definition table.h:36
void SetHasYErr(bool has_yerr)
Definition table.h:81
void set(const Index &i, const double &x, const double &y)
Definition table.h:52
void resize(Index N)
Definition table.cc:38
void Save(std::string filename) const
Definition table.cc:59
int main(int argc, char **argv)
Definition csg_fmatch.cc:39
FileFormatFactory< TrajectoryReader > & TrjReaderFactory()
Eigen::VectorXd linalg_constrained_qrsolve(const Eigen::MatrixXd &A, const Eigen::VectorXd &b, const Eigen::MatrixXd &constr)
solves A*x=b under the constraint B*x = 0
Definition linalg.cc:28
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
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