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