00001 #include "colvarmodule.h"
00002 #include "colvarparse.h"
00003 #include "colvarproxy.h"
00004 #include "colvar.h"
00005 #include "colvarbias.h"
00006 #include "colvarbias_meta.h"
00007 #include "colvarbias_abf.h"
00008
00009
00010 colvarmodule::colvarmodule (char const *config_filename,
00011 colvarproxy *proxy_in)
00012 {
00013
00014 if (proxy == NULL) {
00015 proxy = proxy_in;
00016 parse = new colvarparse();
00017 } else {
00018 cvm::fatal_error ("Error: trying to allocate twice the collective "
00019 "variable module.\n");
00020 }
00021
00022 cvm::log (cvm::line_marker);
00023 cvm::log ("Initializing the collective variables module, version "+
00024 cvm::to_str(COLVARS_VERSION)+".\n");
00025
00026 #if (defined (COLVARS_STANDALONE) || defined (COLVARS_GMXTOOLS))
00027 cvm::log ("Note: this version is not linked to a MD simulation program, "
00028 "and can only be used to analyse data.\n");
00029 #endif
00030
00031
00032
00033 it = it_restart = 0;
00034
00035
00036 config_s.open (config_filename);
00037 if (!config_s)
00038 cvm::fatal_error ("Error: in opening configuration file \""+
00039 std::string (config_filename)+"\".\n");
00040
00041
00042 std::string conf = "";
00043 {
00044 std::string line;
00045 while (colvarparse::getline_nocomments (config_s, line))
00046 conf.append (line+"\n");
00047
00048 config_s.close();
00049 }
00050
00051 parse->get_keyval (conf, "analysis", b_analysis, false);
00052
00053 #if defined (COLVARS_STANDALONE)
00054 parse->get_keyval (conf, "timeStep", dt, 1.0);
00055 #endif
00056
00057 if (cvm::debug())
00058 parse->get_keyval (conf, "debugGradientsStepSize", debug_gradients_step_size, 1.0e-03,
00059 colvarparse::parse_silent);
00060
00061 parse->get_keyval (conf, "eigenvalueCrossingThreshold",
00062 colvarmodule::rotation::crossing_threshold, 1.0e-04,
00063 colvarparse::parse_silent);
00064
00065 parse->get_keyval (conf, "colvarsTrajFrequency", cv_traj_freq, 100);
00066 parse->get_keyval (conf, "colvarsRestartFrequency", restart_out_freq,
00067 proxy->restart_frequency());
00068
00069
00070 bool cv_traj_append;
00071 parse->get_keyval (conf, "colvarsTrajAppend", cv_traj_append, false,
00072 colvarparse::parse_silent);
00073
00074
00075 restart_in_name = proxy->input_prefix().size() ?
00076 std::string (proxy->input_prefix()+".colvars.state") :
00077 std::string ("") ;
00078
00079
00080 restart_out_name = proxy->restart_output_prefix().size() ?
00081 std::string (proxy->restart_output_prefix()+".colvars.state") :
00082 std::string ("");
00083
00084 if (restart_out_name.size())
00085 cvm::log ("The restart output state file will be \""+restart_out_name+"\".\n");
00086
00087 output_prefix = proxy->output_prefix();
00088
00089 std::string const cv_traj_name =
00090 (output_prefix.size() ?
00091 std::string (output_prefix+".colvars.traj") :
00092 std::string ("colvars.traj"));
00093 cvm::log ("The trajectory file will be \""+
00094 cv_traj_name+"\".\n");
00095
00096 cvm::log ("The final output state file will be \""+
00097 (output_prefix.size() ?
00098 std::string (output_prefix+".colvars.state") :
00099 std::string ("colvars.state"))+"\".\n");
00100
00101
00102 init_colvars (conf);
00103
00104
00105 init_biases (conf);
00106
00107 #if defined (COLVARS_STANDALONE)
00108 if (b_analysis) {
00109
00110
00111 cvm::log (cvm::line_marker);
00112
00113 parse->get_keyval (conf, "readTrajectory",
00114 cv_traj_read_name);
00115 parse->get_keyval (conf, "readBegin",
00116 cv_traj_read_begin, 0);
00117 parse->get_keyval (conf, "readEnd",
00118 cv_traj_read_end, 0);
00119 }
00120 #endif
00121
00122
00123 parse->check_keywords (conf, "colvarmodule");
00124 cvm::log (cvm::line_marker);
00125
00126
00127 if (restart_in_name.size()) {
00128
00129 std::ifstream input_is (restart_in_name.c_str());
00130 if (!input_is.good())
00131 fatal_error ("Error: in opening restart file \""+
00132 std::string (restart_in_name)+"\".\n");
00133 else {
00134 cvm::log ("Restarting from file \""+restart_in_name+"\".\n");
00135 read_restart (input_is);
00136 cvm::log (cvm::line_marker);
00137 }
00138 }
00139
00140
00141 #if (!defined (COLVARS_STANDALONE) && !defined (COLVARS_GMXTOOLS))
00142
00143 if ((!output_prefix.size()) && (!restart_out_name.size())) {
00144 cvm::fatal_error ("Error: neither the final output state file or "
00145 "the output restart file could be defined, exiting.\n");
00146 }
00147 #endif
00148
00149
00150 if (cv_traj_freq && !cvm::b_analysis) {
00151
00152 if (cvm::debug())
00153 cvm::log ("Opening output file \""+cv_traj_name+"\".\n");
00154
00155 if (cv_traj_append) {
00156 cv_traj_os.open (cv_traj_name.c_str(), std::ios::app);
00157 } else {
00158 proxy->backup_file (cv_traj_name.c_str());
00159 cv_traj_os.open (cv_traj_name.c_str(), std::ios::out);
00160 }
00161 cv_traj_os.setf (std::ios::scientific, std::ios::floatfield);
00162 }
00163
00164 cvm::log ("Collective variables module initialized.\n");
00165 cvm::log (cvm::line_marker);
00166 }
00167
00168
00169 std::istream & colvarmodule::read_restart (std::istream &is)
00170 {
00171 {
00172
00173 std::string restart_conf;
00174 if (is >> colvarparse::read_block ("configuration", restart_conf)) {
00175 parse->get_keyval (restart_conf, "step",
00176 it_restart, (size_t) 0,
00177 colvarparse::parse_silent);
00178 it = it_restart;
00179 }
00180 is.clear();
00181 }
00182
00183
00184 cvm::increase_depth();
00185 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00186 cvi != colvars.end();
00187 cvi++) {
00188 if ( !((*cvi)->read_restart (is)) )
00189 cvm::fatal_error ("Error: in reading restart configuration for collective variable \""+
00190 (*cvi)->name+"\".\n");
00191 }
00192
00193
00194 for (std::vector<colvarbias *>::iterator bi = biases.begin();
00195 bi != biases.end();
00196 bi++) {
00197 if (!((*bi)->read_restart (is)))
00198 fatal_error ("ERROR: in reading restart configuration for bias \""+
00199 (*bi)->name+"\".\n");
00200 }
00201 cvm::decrease_depth();
00202
00203 return is;
00204 }
00205
00206
00207
00208 void colvarmodule::init_colvars (std::string const &conf)
00209 {
00210 if (cvm::debug())
00211 cvm::log ("Initializing the collective variables.\n");
00212
00213 std::string colvar_conf = "";
00214 size_t pos = 0;
00215 while (parse->key_lookup (conf, "colvar", colvar_conf, pos)) {
00216
00217 if (colvar_conf.size()) {
00218 cvm::log (cvm::line_marker);
00219 cvm::increase_depth();
00220 colvars.push_back (new colvar (colvar_conf));
00221 (colvars.back())->check_keywords (colvar_conf, "colvar");
00222 cvm::decrease_depth();
00223 } else {
00224 cvm::log ("Warning: \"colvar\" keyword found without any configuration.\n");
00225 }
00226 colvar_conf = "";
00227 }
00228
00229
00230 if (!colvars.size())
00231 cvm::fatal_error ("Error: no collective variables defined.\n");
00232
00233 if (colvars.size())
00234 cvm::log (cvm::line_marker);
00235 cvm::log ("Collective variables initialized, "+
00236 cvm::to_str (colvars.size())+
00237 " in total.\n");
00238 }
00239
00240
00241 void colvarmodule::init_biases (std::string const &conf)
00242 {
00243 if (cvm::debug())
00244 cvm::log ("Initializing the collective variables biases.\n");
00245
00246 {
00248 std::string abf_conf = "";
00249 size_t abf_pos = 0;
00250 while (parse->key_lookup (conf, "abf", abf_conf, abf_pos)) {
00251 if (abf_conf.size()) {
00252 cvm::log (cvm::line_marker);
00253 cvm::increase_depth();
00254 biases.push_back (new colvarbias_abf (abf_conf, "abf"));
00255 (biases.back())->check_keywords (abf_conf, "abf");
00256 cvm::decrease_depth();
00257 n_abf_biases++;
00258 } else {
00259 cvm::log ("Warning: \"abf\" keyword found without configuration.\n");
00260 }
00261 abf_conf = "";
00262 }
00263 }
00264
00265 {
00267 std::string harm_conf = "";
00268 size_t harm_pos = 0;
00269 while (parse->key_lookup (conf, "harmonic", harm_conf, harm_pos)) {
00270 if (harm_conf.size()) {
00271 cvm::log (cvm::line_marker);
00272 cvm::increase_depth();
00273 biases.push_back (new colvarbias_harmonic (harm_conf, "harmonic"));
00274 (biases.back())->check_keywords (harm_conf, "harmonic");
00275 cvm::decrease_depth();
00276 n_harm_biases++;
00277 } else {
00278 cvm::log ("Warning: \"harmonic\" keyword found without configuration.\n");
00279 }
00280 harm_conf = "";
00281 }
00282 }
00283
00284 {
00286 std::string histo_conf = "";
00287 size_t histo_pos = 0;
00288 while (parse->key_lookup (conf, "histogram", histo_conf, histo_pos)) {
00289 if (histo_conf.size()) {
00290 cvm::log (cvm::line_marker);
00291 cvm::increase_depth();
00292 biases.push_back (new colvarbias_histogram (histo_conf, "histogram"));
00293 (biases.back())->check_keywords (histo_conf, "histogram");
00294 cvm::decrease_depth();
00295 n_histo_biases++;
00296 } else {
00297 cvm::log ("Warning: \"histogram\" keyword found without configuration.\n");
00298 }
00299 histo_conf = "";
00300 }
00301 }
00302
00303 {
00305 std::string meta_conf = "";
00306 size_t meta_pos = 0;
00307 while (parse->key_lookup (conf, "metadynamics", meta_conf, meta_pos)) {
00308 if (meta_conf.size()) {
00309 cvm::log (cvm::line_marker);
00310 cvm::increase_depth();
00311 biases.push_back (new colvarbias_meta (meta_conf, "metadynamics"));
00312 (biases.back())->check_keywords (meta_conf, "metadynamics");
00313 cvm::decrease_depth();
00314 n_meta_biases++;
00315 } else {
00316 cvm::log ("Warning: \"metadynamics\" keyword found without configuration.\n");
00317 }
00318 meta_conf = "";
00319 }
00320 }
00321
00322 if (biases.size())
00323 cvm::log (cvm::line_marker);
00324 cvm::log ("Collective variables biases initialized, "+
00325 cvm::to_str (biases.size())+" in total.\n");
00326 }
00327
00328
00329 void colvarmodule::calc() {
00330
00331 if (cvm::debug()) {
00332 cvm::log (cvm::line_marker);
00333 cvm::log ("Collective variables module, step no. "+
00334 cvm::to_str (cvm::step_relative())+"\n");
00335 }
00336
00337
00338 if (cvm::debug())
00339 cvm::log ("Calculating collective variables.\n");
00340 cvm::increase_depth();
00341 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00342 cvi != colvars.end();
00343 cvi++) {
00344 (*cvi)->calc();
00345 }
00346 cvm::decrease_depth();
00347
00348
00349
00350 if (cvm::debug() && biases.size())
00351 cvm::log ("Updating collective variable biases.\n");
00352 cvm::increase_depth();
00353 for (std::vector<colvarbias *>::iterator bi = biases.begin();
00354 bi != biases.end();
00355 bi++) {
00356 (*bi)->update();
00357 }
00358 cvm::decrease_depth();
00359
00360
00361 if (cvm::debug() && biases.size())
00362 cvm::log ("Collecting forces from all biases.\n");
00363 cvm::increase_depth();
00364 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00365 cvi != colvars.end();
00366 cvi++) {
00367 (*cvi)->reset_bias_force();
00368 }
00369 for (std::vector<colvarbias *>::iterator bi = biases.begin();
00370 bi != biases.end();
00371 bi++) {
00372 (*bi)->communicate_forces();
00373 }
00374 cvm::decrease_depth();
00375
00376 if (cvm::b_analysis) {
00377
00378 if (cvm::debug() && biases.size())
00379 cvm::log ("Perform runtime analyses.\n");
00380 cvm::increase_depth();
00381 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00382 cvi != colvars.end();
00383 cvi++) {
00384 (*cvi)->analyse();
00385 }
00386 for (std::vector<colvarbias *>::iterator bi = biases.begin();
00387 bi != biases.end();
00388 bi++) {
00389 (*bi)->analyse();
00390 }
00391 cvm::decrease_depth();
00392 }
00393
00394
00395
00396 if (cvm::debug())
00397 cvm::log ("Updating the internal degrees of freedom "
00398 "of colvars (if they have any).\n");
00399 cvm::increase_depth();
00400 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00401 cvi != colvars.end();
00402 cvi++) {
00403 (*cvi)->update();
00404 }
00405 cvm::decrease_depth();
00406
00407
00408
00409 if (cvm::debug())
00410 cvm::log ("Communicating forces from the colvars to the atoms.\n");
00411 cvm::increase_depth();
00412 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00413 cvi != colvars.end();
00414 cvi++) {
00415 if ((*cvi)->tasks[colvar::task_gradients])
00416 (*cvi)->communicate_forces();
00417 }
00418 cvm::decrease_depth();
00419
00420
00421 if (restart_out_freq && !cvm::b_analysis) {
00422 if ( (cvm::step_relative() > 0) &&
00423 ((cvm::step_relative() % restart_out_freq) == 0) ) {
00424 cvm::log ("Writing the current state to the restart file.\n");
00425 proxy->backup_file (restart_out_name.c_str());
00426 restart_out_os.open (restart_out_name.c_str());
00427 restart_out_os.setf (std::ios::scientific, std::ios::floatfield);
00428 if (!write_restart (restart_out_os))
00429 cvm::fatal_error ("Error: in writing restart file.\n");
00430 restart_out_os.close();
00431 }
00432 }
00433
00434
00435 if (cv_traj_freq) {
00436
00437 if (cvm::debug())
00438 cvm::log ("Writing trajectory file.\n");
00439
00440
00441 cvm::increase_depth();
00442 if ((cvm::step_relative() % (cv_traj_freq * 1000)) == 0) {
00443 cv_traj_os << "# " << cvm::wrap_string ("step", cvm::it_width-2)
00444 << " ";
00445 if (cvm::debug())
00446 cv_traj_os.flush();
00447 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00448 cvi != colvars.end();
00449 cvi++) {
00450 (*cvi)->write_traj_label (cv_traj_os);
00451 }
00452 cv_traj_os << "\n";
00453 if (cvm::debug())
00454 cv_traj_os.flush();
00455 }
00456 cvm::decrease_depth();
00457
00458
00459 cvm::increase_depth();
00460 if ((cvm::step_relative() % cv_traj_freq) == 0) {
00461 cv_traj_os << std::setw (cvm::it_width) << it
00462 << " ";
00463 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00464 cvi != colvars.end();
00465 cvi++) {
00466 (*cvi)->write_traj (cv_traj_os);
00467 }
00468 cv_traj_os << "\n";
00469 if (cvm::debug())
00470 cv_traj_os.flush();
00471 }
00472 cvm::decrease_depth();
00473
00474 if (restart_out_freq) {
00475
00476 if ( (cvm::step_relative() > 0) &&
00477 ((cvm::step_relative() % restart_out_freq) == 0) ) {
00478 cvm::log ("Synchronizing trajectory file.\n");
00479 cv_traj_os.flush();
00480 }
00481 }
00482 }
00483 }
00484
00485
00486 void colvarmodule::analyse()
00487 {
00488 if (cvm::debug()) {
00489 cvm::log ("colvarmodule::analyse(), step = "+cvm::to_str (it)+".\n");
00490 }
00491
00492 if (cvm::step_relative() == 0)
00493 cvm::log ("Performing analysis.\n");
00494
00495
00496 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00497 cvi != colvars.end();
00498 cvi++) {
00499 cvm::increase_depth();
00500 (*cvi)->analyse();
00501 cvm::decrease_depth();
00502 }
00503
00504
00505 for (std::vector<colvarbias *>::iterator bi = biases.begin();
00506 bi != biases.end();
00507 bi++) {
00508 cvm::increase_depth();
00509 (*bi)->analyse();
00510 cvm::decrease_depth();
00511 }
00512
00513 }
00514
00515
00516 colvarmodule::~colvarmodule()
00517 {
00518 delete parse;
00519 finalise();
00520 }
00521
00522
00523 void colvarmodule::finalise()
00524 {
00525
00526
00527 #if (!defined (COLVARS_STANDALONE) && !defined (COLVARS_GMXTOOLS))
00528 {
00529
00530
00531 std::string const out_name =
00532 (output_prefix.size() ?
00533 std::string (output_prefix+".colvars.state") :
00534 std::string ("colvars.state"));
00535 cvm::log ("Saving collective variables state to \""+out_name+"\".\n");
00536 proxy->backup_file (out_name.c_str());
00537 std::ofstream out (out_name.c_str());
00538 out.setf (std::ios::scientific, std::ios::floatfield);
00539 this->write_restart (out);
00540 out.close();
00541 }
00542 #endif
00543
00544 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00545 cvi != colvars.end();
00546 cvi++) {
00547 delete *cvi;
00548 }
00549 colvars.clear();
00550
00551 for (std::vector<colvarbias *>::iterator bi = biases.begin();
00552 bi != biases.end();
00553 bi++) {
00554 delete *bi;
00555 }
00556 biases.clear();
00557
00558 if (cv_traj_freq) {
00559 cv_traj_os.close();
00560 }
00561 }
00562
00563
00564
00565 bool colvarmodule::read_traj (char const *traj_filename)
00566 {
00567 std::ifstream static *traj_s = NULL;
00568
00569 if (!traj_s) {
00570 cvm::log ("Opening trajectory file \""+
00571 std::string (traj_filename)+"\".\n");
00572 traj_s = new std::ifstream (traj_filename);
00573 }
00574
00575 do {
00576
00577 std::string line ("");
00578
00579 do {
00580 if (!colvarparse::getline_nocomments (*traj_s, line)) {
00581 cvm::log ("End of file \""+std::string (traj_filename)+
00582 "\" reached.\n");
00583 traj_s->close();
00584 return false;
00585 }
00586 } while (line.find_first_not_of (colvarparse::white_space) == std::string::npos);
00587
00588 std::istringstream is (line);
00589
00590 if (!(is >> it)) return false;
00591
00592 if ( (it < cv_traj_read_begin) ) {
00593
00594 #if defined (COLVARS_STANDALONE)
00595 if ((it % 1000) == 0)
00596 std::cerr << "Skipping trajectory step " << it
00597 << " \r";
00598 #endif
00599
00600 continue;
00601
00602 } else {
00603
00604 #if defined (COLVARS_STANDALONE)
00605 if ((it % 1000) == 0)
00606 std::cerr << "Reading from trajectory, step = " << it
00607 << " \r";
00608 #endif
00609
00610 if ( (cv_traj_read_end > cv_traj_read_begin) &&
00611 (it > cv_traj_read_end) ) {
00612 #if defined (COLVARS_STANDALONE)
00613 std::cerr << "\n";
00614 #endif
00615 cvm::log ("Reached the end of the trajectory, "
00616 "read_end = "+cvm::to_str (cv_traj_read_end)+"\n");
00617 return false;
00618 }
00619
00620 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00621 cvi != colvars.end();
00622 cvi++) {
00623 if (!(*cvi)->read_traj (is)) {
00624 cvm::log ("Error: in reading colvar \""+(*cvi)->name+
00625 "\" from trajectory file \""+
00626 std::string (traj_filename)+"\".\n");
00627 return false;
00628 }
00629 }
00630
00631 break;
00632 }
00633 } while (true);
00634
00635 return true;
00636 }
00637
00638
00639 std::ostream & colvarmodule::write_restart (std::ostream &os)
00640 {
00641 os << "configuration {\n"
00642 << " step " << std::setw (it_width)
00643 << it << "\n"
00644 << " dt " << dt << "\n"
00645 << "}\n\n";
00646
00647 cvm::increase_depth();
00648 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00649 cvi != colvars.end();
00650 cvi++) {
00651 (*cvi)->write_restart (os);
00652 }
00653
00654 for (std::vector<colvarbias *>::iterator bi = biases.begin();
00655 bi != biases.end();
00656 bi++) {
00657 (*bi)->write_restart (os);
00658 }
00659 cvm::decrease_depth();
00660
00661 return os;
00662 }
00663
00664
00665
00666 void cvm::log (std::string const &message)
00667 {
00668 if (depth > 0)
00669 proxy->log ((std::string (2*depth, ' '))+message);
00670 else
00671 proxy->log (message);
00672 }
00673
00674 void cvm::increase_depth()
00675 {
00676 depth++;
00677 }
00678
00679 void cvm::decrease_depth()
00680 {
00681 if (depth) depth--;
00682 }
00683
00684 void cvm::fatal_error (std::string const &message)
00685 {
00686 proxy->fatal_error (message);
00687 }
00688
00689 void cvm::exit (std::string const &message)
00690 {
00691 proxy->exit (message);
00692 }
00693
00694
00695
00696
00697 std::vector<colvar *> colvarmodule::colvars;
00698 std::vector<colvarbias *> colvarmodule::biases;
00699 size_t colvarmodule::n_abf_biases = 0;
00700 size_t colvarmodule::n_harm_biases = 0;
00701 size_t colvarmodule::n_histo_biases = 0;
00702 size_t colvarmodule::n_meta_biases = 0;
00703 colvarproxy *colvarmodule::proxy = NULL;
00704
00705
00706
00707 cvm::real colvarmodule::dt = 1.0;
00708 cvm::real colvarmodule::debug_gradients_step_size = 1.0e-03;
00709 size_t colvarmodule::it = 0;
00710 size_t colvarmodule::it_restart = 0;
00711 size_t colvarmodule::restart_out_freq = 0;
00712 size_t colvarmodule::cv_traj_freq = 0;
00713 size_t colvarmodule::depth = 0;
00714 bool colvarmodule::b_analysis = false;
00715 cvm::real colvarmodule::rotation::crossing_threshold = 1.0E-04;
00716
00717
00718
00719 std::string colvarmodule::output_prefix = "";
00720 std::string colvarmodule::input_prefix = "";
00721 std::string colvarmodule::restart_in_name = "";
00722
00723
00724
00725 size_t const colvarmodule::it_width = 12;
00726 size_t const colvarmodule::cv_prec = 14;
00727 size_t const colvarmodule::cv_width = 21;
00728 size_t const colvarmodule::en_prec = 14;
00729 size_t const colvarmodule::en_width = 21;
00730 std::string const colvarmodule::line_marker =
00731 "----------------------------------------------------------------------\n";
00732
00733
00734
00735
00736
00737
00738 std::ostream & operator << (std::ostream &os, colvarmodule::rvector const &v)
00739 {
00740 std::streamsize const w = os.width();
00741 std::streamsize const p = os.precision();
00742
00743 os.width (2);
00744 os << "( ";
00745 os.width (w); os.precision (p);
00746 os << v.x << " , ";
00747 os.width (w); os.precision (p);
00748 os << v.y << " , ";
00749 os.width (w); os.precision (p);
00750 os << v.z << " )";
00751 return os;
00752 }
00753
00754
00755 std::istream & operator >> (std::istream &is, colvarmodule::rvector &v)
00756 {
00757 size_t const start_pos = is.tellg();
00758 char sep;
00759 if ( !(is >> sep) || !(sep == '(') ||
00760 !(is >> v.x) || !(is >> sep) || !(sep == ',') ||
00761 !(is >> v.y) || !(is >> sep) || !(sep == ',') ||
00762 !(is >> v.z) || !(is >> sep) || !(sep == ')') ) {
00763 is.clear();
00764 is.seekg (start_pos, std::ios::beg);
00765 is.setstate (std::ios::failbit);
00766 return is;
00767 }
00768 return is;
00769 }
00770
00771
00772
00773 std::ostream & operator << (std::ostream &os, colvarmodule::quaternion const &q)
00774 {
00775 std::streamsize const w = os.width();
00776 std::streamsize const p = os.precision();
00777
00778 os.width (2);
00779 os << "( ";
00780 os.width (w); os.precision (p);
00781 os << q.q0 << " , ";
00782 os.width (w); os.precision (p);
00783 os << q.q1 << " , ";
00784 os.width (w); os.precision (p);
00785 os << q.q2 << " , ";
00786 os.width (w); os.precision (p);
00787 os << q.q3 << " )";
00788 return os;
00789 }
00790
00791
00792 std::istream & operator >> (std::istream &is, colvarmodule::quaternion &q)
00793 {
00794 size_t const start_pos = is.tellg();
00795
00796 std::string euler ("");
00797
00798 if ( (is >> euler) && (colvarparse::to_lower_cppstr (euler) ==
00799 std::string ("euler")) ) {
00800
00801
00802
00803 char sep;
00804 cvm::real phi, theta, psi;
00805 if ( !(is >> sep) || !(sep == '(') ||
00806 !(is >> phi) || !(is >> sep) || !(sep == ',') ||
00807 !(is >> theta) || !(is >> sep) || !(sep == ',') ||
00808 !(is >> psi) || !(is >> sep) || !(sep == ')') ) {
00809 is.clear();
00810 is.seekg (start_pos, std::ios::beg);
00811 is.setstate (std::ios::failbit);
00812 return is;
00813 }
00814
00815 q = colvarmodule::quaternion (phi, theta, psi);
00816
00817 } else {
00818
00819
00820
00821 is.seekg (start_pos, std::ios::beg);
00822 char sep;
00823 if ( !(is >> sep) || !(sep == '(') ||
00824 !(is >> q.q0) || !(is >> sep) || !(sep == ',') ||
00825 !(is >> q.q1) || !(is >> sep) || !(sep == ',') ||
00826 !(is >> q.q2) || !(is >> sep) || !(sep == ',') ||
00827 !(is >> q.q3) || !(is >> sep) || !(sep == ')') ) {
00828 is.clear();
00829 is.seekg (start_pos, std::ios::beg);
00830 is.setstate (std::ios::failbit);
00831 return is;
00832 }
00833 }
00834
00835 return is;
00836 }
00837
00838
00839
00840
00841
00842
00843
00844
00845
00846
00847
00848
00849
00850 void colvarmodule::rotation::build_matrix (std::vector<cvm::atom_pos> const &pos1,
00851 std::vector<cvm::atom_pos> const &pos2,
00852 matrix2d<cvm::real, 4, 4> &S)
00853 {
00854 cvm::rmatrix C;
00855
00856
00857 C.reset();
00858 for (size_t i = 0; i < pos1.size(); i++) {
00859 C.xx() += pos1[i].x * pos2[i].x;
00860 C.xy() += pos1[i].x * pos2[i].y;
00861 C.xz() += pos1[i].x * pos2[i].z;
00862 C.yx() += pos1[i].y * pos2[i].x;
00863 C.yy() += pos1[i].y * pos2[i].y;
00864 C.yz() += pos1[i].y * pos2[i].z;
00865 C.zx() += pos1[i].z * pos2[i].x;
00866 C.zy() += pos1[i].z * pos2[i].y;
00867 C.zz() += pos1[i].z * pos2[i].z;
00868 }
00869
00870
00871
00872 S[0][0] = C.xx() + C.yy() + C.zz();
00873 S[1][0] = C.yz() - C.zy();
00874 S[0][1] = S[1][0];
00875 S[2][0] = - C.xz() + C.zx() ;
00876 S[0][2] = S[2][0];
00877 S[3][0] = C.xy() - C.yx();
00878 S[0][3] = S[3][0];
00879 S[1][1] = C.xx() - C.yy() - C.zz();
00880 S[2][1] = C.xy() + C.yx();
00881 S[1][2] = S[2][1];
00882 S[3][1] = C.xz() + C.zx();
00883 S[1][3] = S[3][1];
00884 S[2][2] = - C.xx() + C.yy() - C.zz();
00885 S[3][2] = C.yz() + C.zy();
00886 S[2][3] = S[3][2];
00887 S[3][3] = - C.xx() - C.yy() + C.zz();
00888
00889
00890
00891
00892
00893
00894
00895
00896
00897
00898
00899 }
00900
00901
00902 void colvarmodule::rotation::diagonalize_matrix (matrix2d<cvm::real, 4, 4> &S,
00903 cvm::real S_eigval[4],
00904 matrix2d<cvm::real, 4, 4> &S_eigvec)
00905 {
00906
00907 int jac_nrot = 0;
00908 jacobi (S, 4, S_eigval, S_eigvec, &jac_nrot);
00909 eigsrt (S_eigval, S_eigvec, 4);
00910
00911 transpose (S_eigvec, 4);
00912
00913
00914 for (size_t ie = 0; ie < 4; ie++) {
00915 cvm::real norm2 = 0.0;
00916 for (size_t i = 0; i < 4; i++) norm2 += ::pow (S_eigvec[ie][i], int (2));
00917 cvm::real const norm = ::sqrt (norm2);
00918 for (size_t i = 0; i < 4; i++) S_eigvec[ie][i] /= norm;
00919 }
00920 }
00921
00922
00923
00924
00925 void colvarmodule::rotation::calc_optimal_rotation
00926 (std::vector<cvm::atom_pos> const &pos1,
00927 std::vector<cvm::atom_pos> const &pos2)
00928 {
00929 matrix2d<cvm::real, 4, 4> S;
00930 matrix2d<cvm::real, 4, 4> S_backup;
00931 cvm::real S_eigval[4];
00932 matrix2d<cvm::real, 4, 4> S_eigvec;
00933
00934
00935
00936
00937
00938
00939
00940
00941
00942
00943
00944
00945
00946
00947
00948
00949
00950 build_matrix (pos1, pos2, S);
00951 S_backup = S;
00952
00953 if (cvm::debug()) {
00954 if (b_debug_gradients) {
00955 cvm::log ("S = "+cvm::to_str (cvm::to_str (S_backup), cvm::cv_width, cvm::cv_prec)+"\n");
00956 }
00957 }
00958
00959 diagonalize_matrix (S, S_eigval, S_eigvec);
00960
00961
00962 cvm::real const &L0 = S_eigval[0];
00963 cvm::real const &L1 = S_eigval[1];
00964 cvm::real const &L2 = S_eigval[2];
00965 cvm::real const &L3 = S_eigval[3];
00966 cvm::real const *Q0 = S_eigvec[0];
00967 cvm::real const *Q1 = S_eigvec[1];
00968 cvm::real const *Q2 = S_eigvec[2];
00969 cvm::real const *Q3 = S_eigvec[3];
00970
00971 lambda = L0;
00972 q = cvm::quaternion (Q0);
00973
00974 #if (!defined (COLVARS_STANDALONE) && !defined (COLVARS_GMXTOOLS))
00975 if (q_old.norm2() > 0.0) {
00976 q.match (q_old);
00977 if (q_old.inner (q) < (1.0 - crossing_threshold)) {
00978 cvm::log ("Warning: discontinuous rotation!\n");
00979 }
00980 }
00981 q_old = q;
00982 #endif
00983
00984 if (cvm::debug()) {
00985 if (b_debug_gradients) {
00986 cvm::log ("L0 = "+cvm::to_str (L0, cvm::cv_width, cvm::cv_prec)+
00987 ", Q0 = "+cvm::to_str (cvm::quaternion (Q0), cvm::cv_width, cvm::cv_prec)+
00988 ", Q0*Q0 = "+cvm::to_str (cvm::quaternion (Q0).inner (cvm::quaternion (Q0)), cvm::cv_width, cvm::cv_prec)+
00989 "\n");
00990 cvm::log ("L1 = "+cvm::to_str (L1, cvm::cv_width, cvm::cv_prec)+
00991 ", Q1 = "+cvm::to_str (cvm::quaternion (Q1), cvm::cv_width, cvm::cv_prec)+
00992 ", Q0*Q1 = "+cvm::to_str (cvm::quaternion (Q0).inner (cvm::quaternion (Q1)), cvm::cv_width, cvm::cv_prec)+
00993 "\n");
00994 cvm::log ("L2 = "+cvm::to_str (L2, cvm::cv_width, cvm::cv_prec)+
00995 ", Q2 = "+cvm::to_str (cvm::quaternion (Q2), cvm::cv_width, cvm::cv_prec)+
00996 ", Q0*Q2 = "+cvm::to_str (cvm::quaternion (Q0).inner (cvm::quaternion (Q2)), cvm::cv_width, cvm::cv_prec)+
00997 "\n");
00998 cvm::log ("L3 = "+cvm::to_str (L3, cvm::cv_width, cvm::cv_prec)+
00999 ", Q3 = "+cvm::to_str (cvm::quaternion (Q3), cvm::cv_width, cvm::cv_prec)+
01000 ", Q0*Q3 = "+cvm::to_str (cvm::quaternion (Q0).inner (cvm::quaternion (Q3)), cvm::cv_width, cvm::cv_prec)+
01001 "\n");
01002 }
01003 }
01004
01005
01006
01007
01008 for (size_t ia = 0; ia < dS_1.size(); ia++) {
01009
01010 cvm::real const &a2x = pos2[ia].x;
01011 cvm::real const &a2y = pos2[ia].y;
01012 cvm::real const &a2z = pos2[ia].z;
01013
01014 matrix2d<cvm::rvector, 4, 4> &ds_1 = dS_1[ia];
01015
01016
01017 ds_1.reset();
01018 ds_1[0][0] = cvm::rvector ( a2x, a2y, a2z);
01019 ds_1[1][0] = cvm::rvector ( 0.0, a2z, -a2y);
01020 ds_1[0][1] = ds_1[1][0];
01021 ds_1[2][0] = cvm::rvector (-a2z, 0.0, a2x);
01022 ds_1[0][2] = ds_1[2][0];
01023 ds_1[3][0] = cvm::rvector ( a2y, -a2x, 0.0);
01024 ds_1[0][3] = ds_1[3][0];
01025 ds_1[1][1] = cvm::rvector ( a2x, -a2y, -a2z);
01026 ds_1[2][1] = cvm::rvector ( a2y, a2x, 0.0);
01027 ds_1[1][2] = ds_1[2][1];
01028 ds_1[3][1] = cvm::rvector ( a2z, 0.0, a2x);
01029 ds_1[1][3] = ds_1[3][1];
01030 ds_1[2][2] = cvm::rvector (-a2x, a2y, -a2z);
01031 ds_1[3][2] = cvm::rvector ( 0.0, a2z, a2y);
01032 ds_1[2][3] = ds_1[3][2];
01033 ds_1[3][3] = cvm::rvector (-a2x, -a2y, a2z);
01034
01035 cvm::rvector &dl0_1 = dL0_1[ia];
01036 vector1d<cvm::rvector, 4> &dq0_1 = dQ0_1[ia];
01037
01038
01039
01040
01041
01042 dl0_1.reset();
01043 for (size_t i = 0; i < 4; i++) {
01044 for (size_t j = 0; j < 4; j++) {
01045 dl0_1 += Q0[i] * ds_1[i][j] * Q0[j];
01046 }
01047 }
01048
01049 dq0_1.reset();
01050 for (size_t p = 0; p < 4; p++) {
01051 for (size_t i = 0; i < 4; i++) {
01052 for (size_t j = 0; j < 4; j++) {
01053 dq0_1[p] +=
01054 (Q1[i] * ds_1[i][j] * Q0[j]) / (L0-L1) * Q1[p] +
01055 (Q2[i] * ds_1[i][j] * Q0[j]) / (L0-L2) * Q2[p] +
01056 (Q3[i] * ds_1[i][j] * Q0[j]) / (L0-L3) * Q3[p];
01057 }
01058 }
01059 }
01060 }
01061
01062
01063 for (size_t ia = 0; ia < dS_2.size(); ia++) {
01064
01065 cvm::real const &a1x = pos1[ia].x;
01066 cvm::real const &a1y = pos1[ia].y;
01067 cvm::real const &a1z = pos1[ia].z;
01068
01069 matrix2d<cvm::rvector, 4, 4> &ds_2 = dS_2[ia];
01070
01071 ds_2.reset();
01072 ds_2[0][0] = cvm::rvector ( a1x, a1y, a1z);
01073 ds_2[1][0] = cvm::rvector ( 0.0, -a1z, a1y);
01074 ds_2[0][1] = ds_2[1][0];
01075 ds_2[2][0] = cvm::rvector ( a1z, 0.0, -a1x);
01076 ds_2[0][2] = ds_2[2][0];
01077 ds_2[3][0] = cvm::rvector (-a1y, a1x, 0.0);
01078 ds_2[0][3] = ds_2[3][0];
01079 ds_2[1][1] = cvm::rvector ( a1x, -a1y, -a1z);
01080 ds_2[2][1] = cvm::rvector ( a1y, a1x, 0.0);
01081 ds_2[1][2] = ds_2[2][1];
01082 ds_2[3][1] = cvm::rvector ( a1z, 0.0, a1x);
01083 ds_2[1][3] = ds_2[3][1];
01084 ds_2[2][2] = cvm::rvector (-a1x, a1y, -a1z);
01085 ds_2[3][2] = cvm::rvector ( 0.0, a1z, a1y);
01086 ds_2[2][3] = ds_2[3][2];
01087 ds_2[3][3] = cvm::rvector (-a1x, -a1y, a1z);
01088
01089 cvm::rvector &dl0_2 = dL0_2[ia];
01090 vector1d<cvm::rvector, 4> &dq0_2 = dQ0_2[ia];
01091
01092 dl0_2.reset();
01093 for (size_t i = 0; i < 4; i++) {
01094 for (size_t j = 0; j < 4; j++) {
01095 dl0_2 += Q0[i] * ds_2[i][j] * Q0[j];
01096 }
01097 }
01098
01099 dq0_2.reset();
01100 for (size_t p = 0; p < 4; p++) {
01101 for (size_t i = 0; i < 4; i++) {
01102 for (size_t j = 0; j < 4; j++) {
01103 dq0_2[p] +=
01104 (Q1[i] * ds_2[i][j] * Q0[j]) / (L0-L1) * Q1[p] +
01105 (Q2[i] * ds_2[i][j] * Q0[j]) / (L0-L2) * Q2[p] +
01106 (Q3[i] * ds_2[i][j] * Q0[j]) / (L0-L3) * Q3[p];
01107 }
01108 }
01109 }
01110
01111 if (cvm::debug()) {
01112
01113 if (b_debug_gradients) {
01114
01115 matrix2d<cvm::real, 4, 4> S_new;
01116 cvm::real S_new_eigval[4];
01117 matrix2d<cvm::real, 4, 4> S_new_eigvec;
01118
01119
01120
01121 for (size_t comp = 0; comp < 3; comp++) {
01122
01123 S_new = S_backup;
01124
01125 for (size_t i = 0; i < 4; i++) {
01126 for (size_t j = 0; j < 4; j++) {
01127 S_new[i][j] +=
01128 colvarmodule::debug_gradients_step_size * ds_2[i][j][comp];
01129 }
01130 }
01131
01132
01133
01134 diagonalize_matrix (S_new, S_new_eigval, S_new_eigvec);
01135
01136 cvm::real const &L0_new = S_new_eigval[0];
01137 cvm::real const *Q0_new = S_new_eigvec[0];
01138
01139 cvm::real const DL0 = (dl0_2[comp]) * colvarmodule::debug_gradients_step_size;
01140 cvm::quaternion const q0 (Q0);
01141 cvm::quaternion const DQ0 (dq0_2[0][comp] * colvarmodule::debug_gradients_step_size,
01142 dq0_2[1][comp] * colvarmodule::debug_gradients_step_size,
01143 dq0_2[2][comp] * colvarmodule::debug_gradients_step_size,
01144 dq0_2[3][comp] * colvarmodule::debug_gradients_step_size);
01145
01146 cvm::log ( "|(l_0+dl_0) - l_0^new|/l_0 = "+
01147 cvm::to_str (::fabs (L0+DL0 - L0_new)/L0, cvm::cv_width, cvm::cv_prec)+
01148 ", |(q_0+dq_0) - q_0^new| = "+
01149 cvm::to_str ((Q0+DQ0 - Q0_new).norm(), cvm::cv_width, cvm::cv_prec)+
01150 "\n");
01151 }
01152 }
01153 }
01154 }
01155 }
01156
01157
01158
01159
01160
01161 #define ROTATE(a,i,j,k,l) g=a[i][j];h=a[k][l];a[i][j]=g-s*(h+g*tau); \
01162 a[k][l]=h+s*(g-h*tau);
01163 void jacobi(cvm::real **a, int n, cvm::real d[], cvm::real **v, int *nrot)
01164 {
01165 int j,iq,ip,i;
01166 cvm::real tresh,theta,tau,t,sm,s,h,g,c;
01167
01168 std::vector<cvm::real> b (n, 0.0);
01169 std::vector<cvm::real> z (n, 0.0);
01170
01171 for (ip=0;ip<n;ip++) {
01172 for (iq=0;iq<n;iq++) v[ip][iq]=0.0;
01173 v[ip][ip]=1.0;
01174 }
01175 for (ip=0;ip<n;ip++) {
01176 b[ip]=d[ip]=a[ip][ip];
01177 z[ip]=0.0;
01178 }
01179 *nrot=0;
01180 for (i=0;i<=50;i++) {
01181 sm=0.0;
01182 for (ip=0;ip<n-1;ip++) {
01183 for (iq=ip+1;iq<n;iq++)
01184 sm += fabs(a[ip][iq]);
01185 }
01186 if (sm == 0.0) {
01187 return;
01188 }
01189 if (i < 4)
01190 tresh=0.2*sm/(n*n);
01191 else
01192 tresh=0.0;
01193 for (ip=0;ip<n-1;ip++) {
01194 for (iq=ip+1;iq<n;iq++) {
01195 g=100.0*fabs(a[ip][iq]);
01196 if (i > 4 && (cvm::real)(fabs(d[ip])+g) == (cvm::real)fabs(d[ip])
01197 && (cvm::real)(fabs(d[iq])+g) == (cvm::real)fabs(d[iq]))
01198 a[ip][iq]=0.0;
01199 else if (fabs(a[ip][iq]) > tresh) {
01200 h=d[iq]-d[ip];
01201 if ((cvm::real)(fabs(h)+g) == (cvm::real)fabs(h))
01202 t=(a[ip][iq])/h;
01203 else {
01204 theta=0.5*h/(a[ip][iq]);
01205 t=1.0/(fabs(theta)+sqrt(1.0+theta*theta));
01206 if (theta < 0.0) t = -t;
01207 }
01208 c=1.0/sqrt(1+t*t);
01209 s=t*c;
01210 tau=s/(1.0+c);
01211 h=t*a[ip][iq];
01212 z[ip] -= h;
01213 z[iq] += h;
01214 d[ip] -= h;
01215 d[iq] += h;
01216 a[ip][iq]=0.0;
01217 for (j=0;j<=ip-1;j++) {
01218 ROTATE(a,j,ip,j,iq)
01219 }
01220 for (j=ip+1;j<=iq-1;j++) {
01221 ROTATE(a,ip,j,j,iq)
01222 }
01223 for (j=iq+1;j<n;j++) {
01224 ROTATE(a,ip,j,iq,j)
01225 }
01226 for (j=0;j<n;j++) {
01227 ROTATE(v,j,ip,j,iq)
01228 }
01229 ++(*nrot);
01230 }
01231 }
01232 }
01233 for (ip=0;ip<n;ip++) {
01234 b[ip] += z[ip];
01235 d[ip]=b[ip];
01236 z[ip]=0.0;
01237 }
01238 }
01239 cvm::fatal_error ("Too many iterations in routine jacobi.\n");
01240 }
01241 #undef ROTATE
01242
01243
01244 void eigsrt(cvm::real d[], cvm::real **v, int n)
01245 {
01246 int k,j,i;
01247 cvm::real p;
01248
01249 for (i=0;i<n;i++) {
01250 p=d[k=i];
01251 for (j=i+1;j<n;j++)
01252 if (d[j] >= p) p=d[k=j];
01253 if (k != i) {
01254 d[k]=d[i];
01255 d[i]=p;
01256 for (j=0;j<n;j++) {
01257 p=v[j][i];
01258 v[j][i]=v[j][k];
01259 v[j][k]=p;
01260 }
01261 }
01262 }
01263 }
01264
01265
01266 void transpose(cvm::real **v, int n)
01267 {
01268 cvm::real p;
01269 for (int i=0;i<n;i++) {
01270 for (int j=i+1;j<n;j++) {
01271 p=v[i][j];
01272 v[i][j]=v[j][i];
01273 v[j][i]=p;
01274 }
01275 }
01276 }
01277