00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include <iostream>
00011 #include <iomanip>
00012 #include <sstream>
00013 #include <fstream>
00014 #include <cstring>
00015 #include <vector>
00016 #include <map>
00017
00018 #include "colvarmodule.h"
00019 #include "colvarparse.h"
00020 #include "colvarproxy.h"
00021 #include "colvar.h"
00022 #include "colvarbias.h"
00023 #include "colvarbias_abf.h"
00024 #include "colvarbias_alb.h"
00025 #include "colvarbias_histogram.h"
00026 #include "colvarbias_histogram_reweight_amd.h"
00027 #include "colvarbias_meta.h"
00028 #include "colvarbias_restraint.h"
00029 #include "colvarscript.h"
00030 #include "colvaratoms.h"
00031 #include "colvarcomp.h"
00032
00033
00034
00036 class colvarmodule::usage {
00037
00038 public:
00039
00041 usage();
00042
00044 int cite_feature(std::string const &feature);
00045
00047 int cite_paper(std::string const &paper);
00048
00050 std::string report(int flag);
00051
00052 protected:
00053
00055 std::map<std::string, int> feature_count_;
00056
00058 std::map<std::string, int> paper_count_;
00059
00061 std::map<std::string, std::string> paper_url_;
00062
00064 std::map<std::string, std::string> paper_bibtex_;
00065
00067 std::map<std::string, std::string> feature_paper_map_;
00068
00069 };
00070
00071
00072 colvarmodule::colvarmodule(colvarproxy *proxy_in)
00073 {
00074 depth_s = 0;
00075 log_level_ = 10;
00076
00077 xyz_reader_use_count = 0;
00078
00079 num_biases_types_used_ =
00080 reinterpret_cast<void *>(new std::map<std::string, int>());
00081
00082 restart_version_str.clear();
00083 restart_version_int = 0;
00084
00085 usage_ = new usage();
00086 usage_->cite_feature("Colvars module");
00087
00088 if (proxy != NULL) {
00089
00090
00091 cvm::error("Error: trying to allocate the collective "
00092 "variable module twice.\n", COLVARS_BUG_ERROR);
00093 return;
00094 }
00095
00096 proxy = proxy_in;
00097 parse = new colvarparse();
00098 version_int = proxy->get_version_from_string(COLVARS_VERSION);
00099
00100 cvm::log(cvm::line_marker);
00101 cvm::log("Initializing the collective variables module, version "+
00102 version()+".\n");
00103 cvm::log("Please cite Fiorin et al, Mol Phys 2013:\n"
00104 " https://doi.org/10.1080/00268976.2013.813594\n"
00105 "as well as all other papers listed below for individual features used.\n");
00106
00107 if (proxy->smp_enabled() == COLVARS_OK) {
00108 cvm::log("SMP parallelism is enabled; if needed, use \"smp off\" to override this.\n");
00109 }
00110
00111 #if (__cplusplus >= 201103L)
00112 cvm::log("This version was built with the C++11 standard or higher.\n");
00113 #else
00114 cvm::log("This version was built without the C++11 standard: some features are disabled.\n"
00115 "Please see the following link for details:\n"
00116 " https://colvars.github.io/README-c++11.html\n");
00117 #endif
00118
00119
00120
00121
00122
00123 colvarmodule::it = colvarmodule::it_restart = 0;
00124
00125 use_scripted_forces = false;
00126 scripting_after_biases = false;
00127
00128 colvarmodule::debug_gradients_step_size = 1.0e-07;
00129
00130 colvarmodule::rotation::monitor_crossings = false;
00131 colvarmodule::rotation::crossing_threshold = 1.0e-02;
00132
00133 cv_traj_freq = 100;
00134 restart_out_freq = proxy->default_restart_frequency();
00135
00136 cv_traj_write_labels = true;
00137
00138
00139 proxy->script = new colvarscript(proxy, this);
00140 }
00141
00142
00143 colvarmodule * colvarmodule::main()
00144 {
00145 return proxy ? proxy->colvars : NULL;
00146 }
00147
00148
00149 std::vector<colvar *> *colvarmodule::variables()
00150 {
00151 return &colvars;
00152 }
00153
00154
00155 std::vector<colvar *> *colvarmodule::variables_active()
00156 {
00157 return &colvars_active;
00158 }
00159
00160
00161 std::vector<colvar *> *colvarmodule::variables_active_smp()
00162 {
00163 return &colvars_smp;
00164 }
00165
00166
00167 std::vector<int> *colvarmodule::variables_active_smp_items()
00168 {
00169 return &colvars_smp_items;
00170 }
00171
00172
00173 std::vector<colvarbias *> *colvarmodule::biases_active()
00174 {
00175 return &(biases_active_);
00176 }
00177
00178
00179 size_t colvarmodule::size() const
00180 {
00181 return colvars.size() + biases.size();
00182 }
00183
00184
00185 int colvarmodule::read_config_file(char const *config_filename)
00186 {
00187 cvm::log(cvm::line_marker);
00188 cvm::log("Reading new configuration from file \""+
00189 std::string(config_filename)+"\":\n");
00190
00191
00192 std::istream &config_s = proxy->input_stream(config_filename,
00193 "configuration file/string");
00194 if (!config_s) {
00195 return cvm::error("Error: in opening configuration file \""+
00196 std::string(config_filename)+"\".\n",
00197 COLVARS_FILE_ERROR);
00198 }
00199
00200
00201 std::string conf = "";
00202 std::string line;
00203 while (parse->read_config_line(config_s, line)) {
00204
00205 if (line.find_first_not_of(colvarparse::white_space) != std::string::npos)
00206 conf.append(line+"\n");
00207 }
00208 proxy->close_input_stream(config_filename);
00209
00210 return parse_config(conf);
00211 }
00212
00213
00214 int colvarmodule::read_config_string(std::string const &config_str)
00215 {
00216 cvm::log(cvm::line_marker);
00217 cvm::log("Reading new configuration:\n");
00218 std::istringstream new_config_s(config_str);
00219
00220
00221 std::string conf = "";
00222 std::string line;
00223 while (parse->read_config_line(new_config_s, line)) {
00224
00225 if (line.find_first_not_of(colvarparse::white_space) != std::string::npos)
00226 conf.append(line+"\n");
00227 }
00228
00229 return parse_config(conf);
00230 }
00231
00232
00233 std::istream & colvarmodule::getline(std::istream &is, std::string &line)
00234 {
00235 std::string l;
00236 if (std::getline(is, l)) {
00237 size_t const sz = l.size();
00238 if (sz > 0) {
00239 if (l[sz-1] == '\r' ) {
00240
00241 line = l.substr(0, sz-1);
00242 } else {
00243 line = l;
00244 }
00245 } else {
00246 line.clear();
00247 }
00248 }
00249 return is;
00250 }
00251
00252
00253 int colvarmodule::parse_config(std::string &conf)
00254 {
00255
00256 extra_conf.clear();
00257
00258
00259 if (colvarparse::check_braces(conf, 0) != COLVARS_OK) {
00260 return cvm::error("Error: unmatched curly braces in configuration.\n",
00261 COLVARS_INPUT_ERROR);
00262 }
00263
00264
00265 colvarparse::check_ascii(conf);
00266
00267
00268 if (catch_input_errors(parse_global_params(conf))) {
00269 return get_error();
00270 }
00271
00272
00273 if (catch_input_errors(parse_colvars(conf))) {
00274 return get_error();
00275 }
00276
00277
00278 if (catch_input_errors(parse_biases(conf))) {
00279 return get_error();
00280 }
00281
00282
00283 if (catch_input_errors(parse->check_keywords(conf, "colvarmodule"))) {
00284 return get_error();
00285 }
00286
00287
00288 if (extra_conf.size()) {
00289 catch_input_errors(parse_global_params(extra_conf));
00290 catch_input_errors(parse_colvars(extra_conf));
00291 catch_input_errors(parse_biases(extra_conf));
00292 parse->check_keywords(extra_conf, "colvarmodule");
00293 extra_conf.clear();
00294 if (get_error() != COLVARS_OK) return get_error();
00295 }
00296
00297 cvm::log(cvm::line_marker);
00298 cvm::log("Collective variables module (re)initialized.\n");
00299 cvm::log(cvm::line_marker);
00300
00301 if (source_Tcl_script.size() > 0) {
00302 run_tcl_script(source_Tcl_script);
00303 }
00304
00305 return get_error();
00306 }
00307
00308
00309 std::string const & colvarmodule::get_config() const
00310 {
00311 return parse->get_config();
00312 }
00313
00314
00315 int colvarmodule::append_new_config(std::string const &new_conf)
00316 {
00317 extra_conf += new_conf;
00318 return COLVARS_OK;
00319 }
00320
00321
00322 void colvarmodule::config_changed()
00323 {
00324 cv_traj_write_labels = true;
00325 }
00326
00327
00328 int colvarmodule::parse_global_params(std::string const &conf)
00329 {
00330
00331 parse->get_keyval(conf, "logLevel", log_level_, log_level_,
00332 colvarparse::parse_silent);
00333 {
00334 std::string units;
00335 if (parse->get_keyval(conf, "units", units)) {
00336 units = colvarparse::to_lower_cppstr(units);
00337 int error_code = proxy->set_unit_system(units, (colvars.size() != 0));
00338 if (error_code != COLVARS_OK) {
00339 return error_code;
00340 }
00341 }
00342 }
00343
00344 {
00345 std::string index_file_name;
00346 size_t pos = 0;
00347 while (parse->key_lookup(conf, "indexFile", &index_file_name, &pos)) {
00348 cvm::log("# indexFile = \""+index_file_name+"\"\n");
00349 read_index_file(index_file_name.c_str());
00350 index_file_name.clear();
00351 }
00352 }
00353
00354 if (parse->get_keyval(conf, "smp", proxy->b_smp_active, proxy->b_smp_active)) {
00355 if (proxy->b_smp_active == false) {
00356 cvm::log("SMP parallelism has been disabled.\n");
00357 }
00358 }
00359
00360 bool b_analysis = true;
00361 if (parse->get_keyval(conf, "analysis", b_analysis, true,
00362 colvarparse::parse_silent)) {
00363 cvm::log("Warning: keyword \"analysis\" is deprecated: it is now set "
00364 "to true; individual analyses are performed only if requested.");
00365 }
00366
00367 parse->get_keyval(conf, "debugGradientsStepSize", debug_gradients_step_size,
00368 debug_gradients_step_size,
00369 colvarparse::parse_silent);
00370
00371 parse->get_keyval(conf, "monitorEigenvalueCrossing",
00372 colvarmodule::rotation::monitor_crossings,
00373 colvarmodule::rotation::monitor_crossings,
00374 colvarparse::parse_silent);
00375 parse->get_keyval(conf, "eigenvalueCrossingThreshold",
00376 colvarmodule::rotation::crossing_threshold,
00377 colvarmodule::rotation::crossing_threshold,
00378 colvarparse::parse_silent);
00379
00380 parse->get_keyval(conf, "colvarsTrajFrequency", cv_traj_freq, cv_traj_freq);
00381 parse->get_keyval(conf, "colvarsRestartFrequency",
00382 restart_out_freq, restart_out_freq);
00383
00384 parse->get_keyval(conf, "scriptedColvarForces",
00385 use_scripted_forces, use_scripted_forces);
00386
00387 parse->get_keyval(conf, "scriptingAfterBiases",
00388 scripting_after_biases, scripting_after_biases);
00389
00390 #if defined(COLVARS_TCL)
00391 parse->get_keyval(conf, "sourceTclFile", source_Tcl_script);
00392 #endif
00393
00394 return cvm::get_error();
00395 }
00396
00397
00398 int colvarmodule::run_tcl_script(std::string const &filename) {
00399
00400 int result = COLVARS_OK;
00401
00402 #if defined(COLVARS_TCL)
00403 result = proxy->tcl_run_file(filename);
00404 #endif
00405
00406 return result;
00407 }
00408
00409
00410 int colvarmodule::parse_colvars(std::string const &conf)
00411 {
00412 if (cvm::debug())
00413 cvm::log("Initializing the collective variables.\n");
00414
00415 std::string colvar_conf = "";
00416 size_t pos = 0;
00417 while (parse->key_lookup(conf, "colvar", &colvar_conf, &pos)) {
00418
00419 if (colvar_conf.size()) {
00420 cvm::log(cvm::line_marker);
00421 cvm::increase_depth();
00422 colvars.push_back(new colvar());
00423 if (((colvars.back())->init(colvar_conf) != COLVARS_OK) ||
00424 ((colvars.back())->check_keywords(colvar_conf, "colvar") != COLVARS_OK)) {
00425 cvm::log("Error while constructing colvar number " +
00426 cvm::to_str(colvars.size()) + " : deleting.");
00427 delete colvars.back();
00428 return COLVARS_ERROR;
00429 }
00430 cvm::decrease_depth();
00431 } else {
00432 cvm::error("Error: \"colvar\" keyword found without any configuration.\n", COLVARS_INPUT_ERROR);
00433 return COLVARS_ERROR;
00434 }
00435 cvm::decrease_depth();
00436 colvar_conf = "";
00437 }
00438
00439 if (pos > 0) {
00440
00441 config_changed();
00442 }
00443
00444 if (!colvars.size()) {
00445 cvm::log("Warning: no collective variables defined.\n");
00446 }
00447
00448 if (colvars.size())
00449 cvm::log(cvm::line_marker);
00450 cvm::log("Collective variables initialized, "+
00451 cvm::to_str(colvars.size())+
00452 " in total.\n");
00453
00454 return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
00455 }
00456
00457
00458 bool colvarmodule::check_new_bias(std::string &conf, char const *key)
00459 {
00460 if (cvm::get_error() ||
00461 (biases.back()->check_keywords(conf, key) != COLVARS_OK)) {
00462 cvm::log("Error while constructing bias number " +
00463 cvm::to_str(biases.size()) + " : deleting.\n");
00464 delete biases.back();
00465 return true;
00466 }
00467 return false;
00468 }
00469
00470
00471 template <class bias_type>
00472 int colvarmodule::parse_biases_type(std::string const &conf,
00473 char const *keyword)
00474 {
00475
00476 std::string const &type_keyword = colvarparse::to_lower_cppstr(keyword);
00477
00478
00479
00480 std::map<std::string, int> *num_biases_types_used =
00481 reinterpret_cast<std::map<std::string, int> *>(num_biases_types_used_);
00482 if (num_biases_types_used->count(type_keyword) == 0) {
00483 (*num_biases_types_used)[type_keyword] = 0;
00484 }
00485
00486 std::string bias_conf = "";
00487 size_t conf_saved_pos = 0;
00488 while (parse->key_lookup(conf, keyword, &bias_conf, &conf_saved_pos)) {
00489 if (bias_conf.size()) {
00490 cvm::log(cvm::line_marker);
00491 cvm::increase_depth();
00492 int &bias_count = (*num_biases_types_used)[type_keyword];
00493 biases.push_back(new bias_type(type_keyword.c_str()));
00494 bias_count += 1;
00495 biases.back()->rank = bias_count;
00496 biases.back()->init(bias_conf);
00497 if (cvm::check_new_bias(bias_conf, keyword) != COLVARS_OK) {
00498 return COLVARS_ERROR;
00499 }
00500 cvm::decrease_depth();
00501 } else {
00502 cvm::error("Error: keyword \""+std::string(keyword)+"\" found without configuration.\n",
00503 COLVARS_INPUT_ERROR);
00504 return COLVARS_ERROR;
00505 }
00506 bias_conf = "";
00507 }
00508 if (conf_saved_pos > 0) {
00509
00510 config_changed();
00511 }
00512 return COLVARS_OK;
00513 }
00514
00515
00516 int colvarmodule::parse_biases(std::string const &conf)
00517 {
00518 if (cvm::debug())
00519 cvm::log("Initializing the collective variables biases.\n");
00520
00522 parse_biases_type<colvarbias_abf>(conf, "abf");
00523
00525 parse_biases_type<colvarbias_alb>(conf, "ALB");
00526
00528 parse_biases_type<colvarbias_restraint_harmonic>(conf, "harmonic");
00529
00531 parse_biases_type<colvarbias_restraint_harmonic_walls>(conf, "harmonicWalls");
00532
00534 parse_biases_type<colvarbias_histogram>(conf, "histogram");
00535
00537 parse_biases_type<colvarbias_restraint_histogram>(conf, "histogramRestraint");
00538
00540 parse_biases_type<colvarbias_restraint_linear>(conf, "linear");
00541
00543 parse_biases_type<colvarbias_meta>(conf, "metadynamics");
00544
00546 parse_biases_type<colvarbias_reweightaMD>(conf, "reweightaMD");
00547
00548 if (use_scripted_forces) {
00549 cvm::log(cvm::line_marker);
00550 cvm::increase_depth();
00551 cvm::log("User forces script will be run at each bias update.\n");
00552 cvm::decrease_depth();
00553 }
00554
00555 std::vector<std::string> const time_biases = time_dependent_biases();
00556 if (time_biases.size() > 1) {
00557 cvm::log("WARNING: there are "+cvm::to_str(time_biases.size())+
00558 " time-dependent biases with non-zero force parameters:\n"+
00559 cvm::to_str(time_biases)+"\n"+
00560 "Please ensure that their forces do not counteract each other.\n");
00561 }
00562
00563 if (num_biases() || use_scripted_forces) {
00564 cvm::log(cvm::line_marker);
00565 cvm::log("Collective variables biases initialized, "+
00566 cvm::to_str(num_biases())+" in total.\n");
00567 } else {
00568 if (!use_scripted_forces) {
00569 cvm::log("No collective variables biases were defined.\n");
00570 }
00571 }
00572
00573 return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
00574 }
00575
00576
00577 size_t colvarmodule::num_variables() const
00578 {
00579 return colvars.size();
00580 }
00581
00582
00583 size_t colvarmodule::num_variables_feature(int feature_id) const
00584 {
00585 size_t n = 0;
00586 for (std::vector<colvar *>::const_iterator cvi = colvars.begin();
00587 cvi != colvars.end();
00588 cvi++) {
00589 if ((*cvi)->is_enabled(feature_id)) {
00590 n++;
00591 }
00592 }
00593 return n;
00594 }
00595
00596
00597 size_t colvarmodule::num_biases() const
00598 {
00599 return biases.size();
00600 }
00601
00602
00603 size_t colvarmodule::num_biases_feature(int feature_id) const
00604 {
00605 size_t n = 0;
00606 for (std::vector<colvarbias *>::const_iterator bi = biases.begin();
00607 bi != biases.end();
00608 bi++) {
00609 if ((*bi)->is_enabled(feature_id)) {
00610 n++;
00611 }
00612 }
00613 return n;
00614 }
00615
00616
00617 size_t colvarmodule::num_biases_type(std::string const &type) const
00618 {
00619 size_t n = 0;
00620 for (std::vector<colvarbias *>::const_iterator bi = biases.begin();
00621 bi != biases.end();
00622 bi++) {
00623 if ((*bi)->bias_type == type) {
00624 n++;
00625 }
00626 }
00627 return n;
00628 }
00629
00630
00631 std::vector<std::string> const colvarmodule::time_dependent_biases() const
00632 {
00633 size_t i;
00634 std::vector<std::string> biases_names;
00635 for (i = 0; i < num_biases(); i++) {
00636 if (biases[i]->is_enabled(colvardeps::f_cvb_apply_force) &&
00637 biases[i]->is_enabled(colvardeps::f_cvb_active) &&
00638 (biases[i]->is_enabled(colvardeps::f_cvb_history_dependent) ||
00639 biases[i]->is_enabled(colvardeps::f_cvb_time_dependent))) {
00640 biases_names.push_back(biases[i]->name);
00641 }
00642 }
00643 return biases_names;
00644 }
00645
00646
00647 int colvarmodule::catch_input_errors(int result)
00648 {
00649 if (result != COLVARS_OK || get_error()) {
00650 set_error_bits(result);
00651 set_error_bits(COLVARS_INPUT_ERROR);
00652 parse->clear();
00653 return get_error();
00654 }
00655 return COLVARS_OK;
00656 }
00657
00658
00659 colvarbias * colvarmodule::bias_by_name(std::string const &name)
00660 {
00661 colvarmodule *cv = cvm::main();
00662 for (std::vector<colvarbias *>::iterator bi = cv->biases.begin();
00663 bi != cv->biases.end();
00664 bi++) {
00665 if ((*bi)->name == name) {
00666 return (*bi);
00667 }
00668 }
00669 return NULL;
00670 }
00671
00672
00673 colvar *colvarmodule::colvar_by_name(std::string const &name)
00674 {
00675 colvarmodule *cv = cvm::main();
00676 for (std::vector<colvar *>::iterator cvi = cv->colvars.begin();
00677 cvi != cv->colvars.end();
00678 cvi++) {
00679 if ((*cvi)->name == name) {
00680 return (*cvi);
00681 }
00682 }
00683 return NULL;
00684 }
00685
00686
00687 cvm::atom_group *colvarmodule::atom_group_by_name(std::string const &name)
00688 {
00689 colvarmodule *cv = cvm::main();
00690 for (std::vector<cvm::atom_group *>::iterator agi = cv->named_atom_groups.begin();
00691 agi != cv->named_atom_groups.end();
00692 agi++) {
00693 if ((*agi)->name == name) {
00694 return (*agi);
00695 }
00696 }
00697 return NULL;
00698 }
00699
00700
00701 void colvarmodule::register_named_atom_group(atom_group *ag) {
00702 named_atom_groups.push_back(ag);
00703 }
00704
00705
00706 void colvarmodule::unregister_named_atom_group(cvm::atom_group *ag)
00707 {
00708 for (std::vector<cvm::atom_group *>::iterator agi = named_atom_groups.begin();
00709 agi != named_atom_groups.end();
00710 agi++) {
00711 if (*agi == ag) {
00712 named_atom_groups.erase(agi);
00713 break;
00714 }
00715 }
00716 }
00717
00718
00719 int colvarmodule::change_configuration(std::string const &bias_name,
00720 std::string const &conf)
00721 {
00722
00723
00724 cvm::increase_depth();
00725 colvarbias *b;
00726 b = bias_by_name(bias_name);
00727 if (b == NULL) {
00728 cvm::error("Error: bias not found: " + bias_name);
00729 return COLVARS_ERROR;
00730 }
00731 b->change_configuration(conf);
00732 cvm::decrease_depth();
00733 return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
00734 }
00735
00736
00737 std::string colvarmodule::read_colvar(std::string const &name)
00738 {
00739 cvm::increase_depth();
00740 colvar *c;
00741 std::stringstream ss;
00742 c = colvar_by_name(name);
00743 if (c == NULL) {
00744 cvm::error("Error: colvar not found: " + name);
00745 return std::string();
00746 }
00747 ss << c->value();
00748 cvm::decrease_depth();
00749 return ss.str();
00750 }
00751
00752
00753 cvm::real colvarmodule::energy_difference(std::string const &bias_name,
00754 std::string const &conf)
00755 {
00756 cvm::increase_depth();
00757 colvarbias *b;
00758 cvm::real energy_diff = 0.;
00759 b = bias_by_name(bias_name);
00760 if (b == NULL) {
00761 cvm::error("Error: bias not found: " + bias_name);
00762 return 0.;
00763 }
00764 energy_diff = b->energy_difference(conf);
00765 cvm::decrease_depth();
00766 return energy_diff;
00767 }
00768
00769
00770 int colvarmodule::calc()
00771 {
00772 int error_code = COLVARS_OK;
00773
00774 if (cvm::debug()) {
00775 cvm::log(cvm::line_marker);
00776 cvm::log("Collective variables module, step no. "+
00777 cvm::to_str(cvm::step_absolute())+"\n");
00778 }
00779
00780 error_code |= calc_colvars();
00781 error_code |= calc_biases();
00782 error_code |= update_colvar_forces();
00783
00784 error_code |= analyze();
00785
00786
00787 if (cv_traj_freq && cv_traj_name.size()) {
00788 error_code |= write_traj_files();
00789 }
00790
00791
00792 if (restart_out_freq && (cvm::step_relative() > 0) &&
00793 ((cvm::step_absolute() % restart_out_freq) == 0) ) {
00794
00795 if (restart_out_name.size()) {
00796
00797 error_code |= write_restart_file(restart_out_name);
00798 } else {
00799 error_code |= write_restart_file(output_prefix()+".colvars.state");
00800 }
00801
00802 cvm::increase_depth();
00803 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00804 cvi != colvars.end();
00805 cvi++) {
00806
00807 error_code |= (*cvi)->write_output_files();
00808 }
00809 for (std::vector<colvarbias *>::iterator bi = biases.begin();
00810 bi != biases.end();
00811 bi++) {
00812 error_code |= (*bi)->write_state_to_replicas();
00813 }
00814 cvm::decrease_depth();
00815 }
00816
00817
00818 cvm::increase_depth();
00819 for (std::vector<colvarbias *>::iterator bi = biases.begin();
00820 bi != biases.end();
00821 bi++) {
00822 if ((*bi)->output_freq > 0) {
00823 if ((cvm::step_relative() > 0) &&
00824 ((cvm::step_absolute() % (*bi)->output_freq) == 0) ) {
00825 error_code |= (*bi)->write_output_files();
00826 }
00827 }
00828 }
00829 cvm::decrease_depth();
00830
00831 error_code |= end_of_step();
00832
00833
00834 error_code |= proxy->end_of_step();
00835
00836 return error_code;
00837 }
00838
00839
00840 int colvarmodule::calc_colvars()
00841 {
00842 if (cvm::debug())
00843 cvm::log("Calculating collective variables.\n");
00844
00845
00846
00847
00848 std::vector<colvarbias *>::iterator bi;
00849 for (bi = biases.begin(); bi != biases.end(); bi++) {
00850 int const tsf = (*bi)->get_time_step_factor();
00851 if (tsf > 1) {
00852 if (step_absolute() % tsf == 0) {
00853 (*bi)->enable(colvardeps::f_cvb_awake);
00854 } else {
00855 (*bi)->disable(colvardeps::f_cvb_awake);
00856 }
00857 }
00858 }
00859
00860 int error_code = COLVARS_OK;
00861 std::vector<colvar *>::iterator cvi;
00862
00863
00864 variables_active()->clear();
00865 variables_active()->reserve(variables()->size());
00866 for (cvi = variables()->begin(); cvi != variables()->end(); cvi++) {
00867
00868 int tsf = (*cvi)->get_time_step_factor();
00869 if (tsf > 1) {
00870 if (step_absolute() % tsf == 0) {
00871 (*cvi)->enable(colvardeps::f_cv_awake);
00872 } else {
00873 (*cvi)->disable(colvardeps::f_cv_awake);
00874 }
00875 }
00876
00877 if ((*cvi)->is_enabled()) {
00878 variables_active()->push_back(*cvi);
00879 }
00880 }
00881
00882
00883 if (proxy->smp_enabled() == COLVARS_OK) {
00884
00885
00886
00887 variables_active_smp()->clear();
00888 variables_active_smp_items()->clear();
00889
00890 variables_active_smp()->reserve(variables_active()->size());
00891 variables_active_smp_items()->reserve(variables_active()->size());
00892
00893
00894 cvm::increase_depth();
00895 for (cvi = variables_active()->begin(); cvi != variables_active()->end(); cvi++) {
00896
00897 error_code |= (*cvi)->update_cvc_flags();
00898
00899 size_t num_items = (*cvi)->num_active_cvcs();
00900 variables_active_smp()->reserve(variables_active_smp()->size() + num_items);
00901 variables_active_smp_items()->reserve(variables_active_smp_items()->size() + num_items);
00902 for (size_t icvc = 0; icvc < num_items; icvc++) {
00903 variables_active_smp()->push_back(*cvi);
00904 variables_active_smp_items()->push_back(icvc);
00905 }
00906 }
00907 cvm::decrease_depth();
00908
00909
00910 error_code |= proxy->smp_colvars_loop();
00911
00912 cvm::increase_depth();
00913 for (cvi = variables_active()->begin(); cvi != variables_active()->end(); cvi++) {
00914 error_code |= (*cvi)->collect_cvc_data();
00915 }
00916 cvm::decrease_depth();
00917
00918 } else {
00919
00920
00921 cvm::increase_depth();
00922 for (cvi = variables_active()->begin(); cvi != variables_active()->end(); cvi++) {
00923 error_code |= (*cvi)->calc();
00924 if (cvm::get_error()) {
00925 return COLVARS_ERROR;
00926 }
00927 }
00928 cvm::decrease_depth();
00929 }
00930
00931 error_code |= cvm::get_error();
00932 return error_code;
00933 }
00934
00935
00936 int colvarmodule::calc_biases()
00937 {
00938
00939
00940 if (cvm::debug() && num_biases())
00941 cvm::log("Updating collective variable biases.\n");
00942
00943
00944 for (std::vector<colvar *>::iterator cvi = colvars.begin();
00945 cvi != colvars.end(); cvi++) {
00946 (*cvi)->reset_bias_force();
00947 }
00948
00949 std::vector<colvarbias *>::iterator bi;
00950 int error_code = COLVARS_OK;
00951
00952
00953 total_bias_energy = 0.0;
00954
00955
00956
00957 biases_active()->clear();
00958 biases_active()->reserve(biases.size());
00959 for (bi = biases.begin(); bi != biases.end(); bi++) {
00960 if ((*bi)->is_enabled()) {
00961 biases_active()->push_back(*bi);
00962 }
00963 }
00964
00965
00966 if (proxy->smp_enabled() == COLVARS_OK) {
00967
00968 if (use_scripted_forces && !scripting_after_biases) {
00969
00970 error_code |= proxy->smp_biases_script_loop();
00971 } else {
00972
00973 error_code |= proxy->smp_biases_loop();
00974 }
00975
00976 } else {
00977
00978 if (use_scripted_forces && !scripting_after_biases) {
00979 error_code |= calc_scripted_forces();
00980 }
00981
00982 cvm::increase_depth();
00983 for (bi = biases_active()->begin(); bi != biases_active()->end(); bi++) {
00984 error_code |= (*bi)->update();
00985 if (cvm::get_error()) {
00986 return error_code;
00987 }
00988 }
00989 cvm::decrease_depth();
00990 }
00991
00992 for (bi = biases_active()->begin(); bi != biases_active()->end(); bi++) {
00993 total_bias_energy += (*bi)->get_energy();
00994 }
00995
00996 return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
00997 }
00998
00999
01000 int colvarmodule::update_colvar_forces()
01001 {
01002 int error_code = COLVARS_OK;
01003
01004 std::vector<colvar *>::iterator cvi;
01005 std::vector<colvarbias *>::iterator bi;
01006
01007
01008 if (cvm::debug() && num_biases())
01009 cvm::log("Collecting forces from all biases.\n");
01010 cvm::increase_depth();
01011 for (bi = biases_active()->begin(); bi != biases_active()->end(); bi++) {
01012 error_code |= (*bi)->communicate_forces();
01013 }
01014 cvm::decrease_depth();
01015
01016 if (use_scripted_forces && scripting_after_biases) {
01017 error_code |= calc_scripted_forces();
01018 }
01019
01020
01021 if (cvm::debug())
01022 cvm::log("Adding total bias energy: " + cvm::to_str(total_bias_energy) + "\n");
01023 proxy->add_energy(total_bias_energy);
01024
01025 cvm::real total_colvar_energy = 0.0;
01026
01027
01028
01029 if (cvm::debug())
01030 cvm::log("Updating the internal degrees of freedom "
01031 "of colvars (if they have any).\n");
01032 cvm::increase_depth();
01033 for (cvi = variables()->begin(); cvi != variables()->end(); cvi++) {
01034
01035 total_colvar_energy += (*cvi)->update_forces_energy();
01036 }
01037 cvm::decrease_depth();
01038 if (cvm::debug())
01039 cvm::log("Adding total colvar energy: " + cvm::to_str(total_colvar_energy) + "\n");
01040 proxy->add_energy(total_colvar_energy);
01041
01042
01043
01044 if (cvm::debug())
01045 cvm::log("Communicating forces from the colvars to the atoms.\n");
01046 cvm::increase_depth();
01047 for (cvi = variables_active()->begin(); cvi != variables_active()->end(); cvi++) {
01048 if ((*cvi)->is_enabled(colvardeps::f_cv_gradient)) {
01049 (*cvi)->communicate_forces();
01050 if (cvm::get_error()) {
01051 return COLVARS_ERROR;
01052 }
01053 }
01054 }
01055 cvm::decrease_depth();
01056
01057 return error_code;
01058 }
01059
01060
01061 int colvarmodule::calc_scripted_forces()
01062 {
01063
01064
01065 int res;
01066 res = proxy->run_force_callback();
01067 if (res == COLVARS_NOT_IMPLEMENTED) {
01068 cvm::error("Colvar forces scripts are not implemented.");
01069 return COLVARS_NOT_IMPLEMENTED;
01070 }
01071 if (res != COLVARS_OK) {
01072 cvm::error("Error running user colvar forces script");
01073 return COLVARS_ERROR;
01074 }
01075 return COLVARS_OK;
01076 }
01077
01078
01079 int colvarmodule::write_restart_file(std::string const &out_name)
01080 {
01081 cvm::log("Saving collective variables state to \""+out_name+"\".\n");
01082 std::ostream &restart_out_os = proxy->output_stream(out_name, "state file");
01083 if (!restart_out_os) return COLVARS_FILE_ERROR;
01084 if (!write_restart(restart_out_os)) {
01085 return cvm::error("Error: in writing restart file.\n", COLVARS_FILE_ERROR);
01086 }
01087 proxy->close_output_stream(out_name);
01088
01089
01090
01091 return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
01092 }
01093
01094
01095 int colvarmodule::write_restart_string(std::string &output)
01096 {
01097 cvm::log("Saving state to output buffer.\n");
01098 std::ostringstream os;
01099 if (!write_restart(os)) {
01100 return cvm::error("Error: in writing restart to buffer.\n", COLVARS_FILE_ERROR);
01101 }
01102 output = os.str();
01103 return COLVARS_OK;
01104 }
01105
01106
01107 int colvarmodule::write_traj_files()
01108 {
01109 int error_code = COLVARS_OK;
01110
01111 if (cvm::debug()) {
01112 cvm::log("colvarmodule::write_traj_files()\n");
01113 }
01114
01115 std::ostream &cv_traj_os = proxy->output_stream(cv_traj_name,
01116 "colvars trajectory");
01117
01118 if (!cv_traj_os) {
01119 return COLVARS_FILE_ERROR;
01120 }
01121
01122
01123 if ( (cvm::step_relative() == 0) || cv_traj_write_labels ||
01124 ((cvm::step_absolute() % (cv_traj_freq * 1000)) == 0) ) {
01125 error_code |=
01126 write_traj_label(cv_traj_os) ? COLVARS_OK : COLVARS_FILE_ERROR;
01127 cv_traj_write_labels = false;
01128 }
01129
01130 if (cvm::debug()) {
01131 proxy->flush_output_stream(cv_traj_name);
01132 }
01133
01134 if ((cvm::step_absolute() % cv_traj_freq) == 0) {
01135 error_code |= write_traj(cv_traj_os) ? COLVARS_OK : COLVARS_FILE_ERROR;
01136 }
01137
01138 if (cvm::debug()) {
01139 proxy->flush_output_stream(cv_traj_name);
01140 }
01141
01142 if (restart_out_freq && ((cvm::step_absolute() % restart_out_freq) == 0)) {
01143 cvm::log("Synchronizing (emptying the buffer of) trajectory file \""+
01144 cv_traj_name+"\".\n");
01145 error_code |= proxy->flush_output_stream(cv_traj_name);
01146 }
01147
01148 return error_code;
01149 }
01150
01151
01152 int colvarmodule::analyze()
01153 {
01154 if (cvm::debug()) {
01155 cvm::log("colvarmodule::analyze(), step = "+cvm::to_str(it)+".\n");
01156 }
01157
01158
01159 for (std::vector<colvar *>::iterator cvi = variables_active()->begin();
01160 cvi != variables_active()->end();
01161 cvi++) {
01162 cvm::increase_depth();
01163 (*cvi)->analyze();
01164 cvm::decrease_depth();
01165 }
01166
01167
01168 for (std::vector<colvarbias *>::iterator bi = biases.begin();
01169 bi != biases.end();
01170 bi++) {
01171 cvm::increase_depth();
01172 (*bi)->analyze();
01173 cvm::decrease_depth();
01174 }
01175
01176 return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
01177 }
01178
01179
01180 int colvarmodule::end_of_step()
01181 {
01182 if (cvm::debug()) {
01183 cvm::log("colvarmodule::end_of_step(), step = "+cvm::to_str(it)+".\n");
01184 }
01185
01186 for (std::vector<colvar *>::iterator cvi = variables_active()->begin();
01187 cvi != variables_active()->end();
01188 cvi++) {
01189 cvm::increase_depth();
01190 (*cvi)->end_of_step();
01191 cvm::decrease_depth();
01192 }
01193
01194
01195 for (std::vector<colvarbias *>::iterator bi = biases.begin();
01196 bi != biases.end();
01197 bi++) {
01198 cvm::increase_depth();
01199 (*bi)->end_of_step();
01200 cvm::decrease_depth();
01201 }
01202
01203 return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
01204 }
01205
01206
01207 int colvarmodule::update_engine_parameters()
01208 {
01209 if (this->size() == 0) return cvm::get_error();
01210 for (std::vector<colvar *>::iterator cvi = variables()->begin();
01211 cvi != variables()->end(); cvi++) {
01212 (*cvi)->setup();
01213 }
01214 return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
01215 }
01216
01217
01218 colvarmodule::~colvarmodule()
01219 {
01220 if ((proxy->smp_thread_id() < 0) ||
01221 (proxy->smp_thread_id() == 0)) {
01222
01223 reset();
01224
01225
01226 colvarbias::delete_features();
01227 colvar::delete_features();
01228 colvar::cvc::delete_features();
01229 atom_group::delete_features();
01230
01231 delete
01232 reinterpret_cast<std::map<std::string, int> *>(num_biases_types_used_);
01233 num_biases_types_used_ = NULL;
01234
01235 delete parse;
01236 parse = NULL;
01237
01238 delete usage_;
01239 usage_ = NULL;
01240
01241
01242 proxy = NULL;
01243 }
01244 }
01245
01246
01247 int colvarmodule::reset()
01248 {
01249 cvm::log("Resetting the Collective Variables module.\n");
01250
01251 parse->clear();
01252
01253
01254 for (std::vector<colvarbias *>::reverse_iterator bi = biases.rbegin();
01255 bi != biases.rend();
01256 bi++) {
01257 delete *bi;
01258 }
01259 biases.clear();
01260 biases_active_.clear();
01261
01262
01263 reinterpret_cast<std::map<std::string, int> *>(num_biases_types_used_)->clear();
01264
01265
01266 for (std::vector<colvar *>::reverse_iterator cvi = colvars.rbegin();
01267 cvi != colvars.rend();
01268 cvi++) {
01269 delete *cvi;
01270 }
01271 colvars.clear();
01272
01273 reset_index_groups();
01274
01275 proxy->flush_output_streams();
01276 proxy->reset();
01277
01278 return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
01279 }
01280
01281
01282 int colvarmodule::setup_input()
01283 {
01284 if (proxy->input_prefix().size()) {
01285
01286 std::string restart_in_name(proxy->input_prefix()+
01287 std::string(".colvars.state"));
01288 std::istream *input_is = &(proxy->input_stream(restart_in_name,
01289 "restart file/channel",
01290 false));
01291 if (!*input_is) {
01292
01293 restart_in_name = proxy->input_prefix();
01294 input_is = &(proxy->input_stream(restart_in_name,
01295 "restart file/channel"));
01296 if (!*input_is) {
01297 return COLVARS_FILE_ERROR;
01298 }
01299 }
01300
01301
01302
01303 proxy->input_prefix().clear();
01304
01305 cvm::log(cvm::line_marker);
01306 cvm::log("Loading state from file \""+restart_in_name+"\".\n");
01307 read_restart(*input_is);
01308 cvm::log(cvm::line_marker);
01309
01310 proxy->close_input_stream(restart_in_name);
01311
01312 return cvm::get_error();
01313 }
01314
01315
01316 if (proxy->input_buffer() != NULL) {
01317
01318 char const *buffer = proxy->input_buffer();
01319 size_t const buffer_size = strlen(proxy->input_buffer());
01320
01321 proxy->input_buffer() = NULL;
01322 if (buffer_size > 0) {
01323 std::istringstream input_is;
01324
01325
01326 input_is.rdbuf()->pubsetbuf(const_cast<char *>(buffer), buffer_size);
01327 cvm::log(cvm::line_marker);
01328 cvm::log("Loading state from input buffer.\n");
01329 read_restart(input_is);
01330 cvm::log(cvm::line_marker);
01331 return cvm::get_error();
01332 }
01333 }
01334
01335 return COLVARS_OK;
01336 }
01337
01338
01339 int colvarmodule::setup_output()
01340 {
01341 int error_code = COLVARS_OK;
01342
01343
01344 restart_out_name = proxy->restart_output_prefix().size() ?
01345 std::string(proxy->restart_output_prefix()+".colvars.state") :
01346 std::string("");
01347
01348 if (restart_out_name.size()) {
01349 cvm::log("The restart output state file will be \""+
01350 restart_out_name+"\".\n");
01351 }
01352
01353 output_prefix() = proxy->output_prefix();
01354 if (output_prefix().size()) {
01355 cvm::log("The final output state file will be \""+
01356 (output_prefix().size() ?
01357 std::string(output_prefix()+".colvars.state") :
01358 std::string("colvars.state"))+"\".\n");
01359
01360 }
01361
01362 cv_traj_name =
01363 (output_prefix().size() ?
01364 std::string(output_prefix()+".colvars.traj") :
01365 std::string(""));
01366
01367 for (std::vector<colvarbias *>::iterator bi = biases.begin();
01368 bi != biases.end();
01369 bi++) {
01370 error_code |= (*bi)->setup_output();
01371 }
01372
01373 if (error_code != COLVARS_OK || cvm::get_error()) {
01374 set_error_bits(COLVARS_FILE_ERROR);
01375 }
01376
01377 return cvm::get_error();
01378 }
01379
01380
01381 std::string colvarmodule::state_file_prefix(char const *filename)
01382 {
01383 std::string const filename_str(filename);
01384 std::string const prefix =
01385 filename_str.substr(0, filename_str.find(".colvars.state"));
01386 if (prefix.size() == 0) {
01387 cvm::error("Error: invalid filename/prefix value \""+filename_str+"\".",
01388 COLVARS_INPUT_ERROR);
01389 }
01390 return prefix;
01391 }
01392
01393
01394
01395 std::istream & colvarmodule::read_restart(std::istream &is)
01396 {
01397 bool warn_total_forces = false;
01398
01399 {
01400
01401 std::string restart_conf;
01402 if (is >> colvarparse::read_block("configuration", &restart_conf)) {
01403
01404 parse->get_keyval(restart_conf, "step",
01405 it_restart, static_cast<step_number>(0),
01406 colvarparse::parse_restart);
01407 it = it_restart;
01408
01409 restart_version_str.clear();
01410 restart_version_int = 0;
01411 parse->get_keyval(restart_conf, "version",
01412 restart_version_str, std::string(""),
01413 colvarparse::parse_restart);
01414 if (restart_version_str.size()) {
01415
01416 restart_version_int =
01417 proxy->get_version_from_string(restart_version_str.c_str());
01418 }
01419
01420 if (restart_version() != version()) {
01421 cvm::log("This state file was generated with version "+
01422 restart_version()+"\n");
01423 }
01424
01425 if (restart_version_number() < 20160810) {
01426
01427 if (proxy->total_forces_enabled()) {
01428 warn_total_forces = true;
01429 }
01430 }
01431
01432 std::string units_restart;
01433 if (parse->get_keyval(restart_conf, "units",
01434 units_restart, std::string(""),
01435 colvarparse::parse_restart)) {
01436 units_restart = colvarparse::to_lower_cppstr(units_restart);
01437 if ((proxy->units.size() > 0) && (units_restart != proxy->units)) {
01438 cvm::error("Error: the state file has units \""+units_restart+
01439 "\", but the current unit system is \""+proxy->units+
01440 "\".\n", COLVARS_INPUT_ERROR);
01441 }
01442 }
01443
01444 }
01445 is.clear();
01446 parse->clear_keyword_registry();
01447 }
01448
01449 print_total_forces_errning(warn_total_forces);
01450
01451 read_objects_state(is);
01452
01453 return is;
01454 }
01455
01456
01457
01458 std::istream & colvarmodule::read_objects_state(std::istream &is)
01459 {
01460 std::streampos pos = 0;
01461 std::string word;
01462
01463 while (is.good()) {
01464 pos = is.tellg();
01465 word.clear();
01466 is >> word;
01467
01468 if (word.size()) {
01469
01470 is.seekg(pos, std::ios::beg);
01471
01472 if (word == "colvar") {
01473
01474 cvm::increase_depth();
01475 for (std::vector<colvar *>::iterator cvi = colvars.begin();
01476 cvi != colvars.end();
01477 cvi++) {
01478 if ( !((*cvi)->read_state(is)) ) {
01479
01480
01481 cvm::error("Error: in reading restart configuration for "
01482 "collective variable \""+(*cvi)->name+"\".\n",
01483 COLVARS_INPUT_ERROR);
01484 }
01485 if (is.tellg() > pos) break;
01486 }
01487 cvm::decrease_depth();
01488
01489 } else {
01490
01491 cvm::increase_depth();
01492 for (std::vector<colvarbias *>::iterator bi = biases.begin();
01493 bi != biases.end();
01494 bi++) {
01495 if (((*bi)->state_keyword != word) && (*bi)->bias_type != word) {
01496
01497
01498 continue;
01499 }
01500 if (!((*bi)->read_state(is))) {
01501
01502 cvm::error("Error: in reading restart configuration for bias \""+
01503 (*bi)->name+"\".\n",
01504 COLVARS_INPUT_ERROR);
01505 }
01506 if (is.tellg() > pos) break;
01507 }
01508 cvm::decrease_depth();
01509 }
01510 }
01511
01512 if (is.tellg() == pos) {
01513
01514
01515 is >> colvarparse::read_block(word, NULL);
01516 }
01517
01518 if (!is) break;
01519 }
01520
01521 return is;
01522 }
01523
01524
01525 int colvarmodule::print_total_forces_errning(bool warn_total_forces)
01526 {
01527 if (warn_total_forces) {
01528 cvm::log(cvm::line_marker);
01529 cvm::log("WARNING: The definition of system forces has changed. Please see:\n");
01530 cvm::log(" https://colvars.github.io/README-totalforce.html\n");
01531
01532 output_prefix() = proxy->input_prefix();
01533 cvm::log("All output files will now be saved with the prefix \""+output_prefix()+".tmp.*\".\n");
01534 cvm::log("Please review the important warning above. After that, you may rename:\n\
01535 \""+output_prefix()+".tmp.colvars.state\"\n\
01536 to:\n\
01537 \""+proxy->input_prefix()+".colvars.state\"\n\
01538 and load it to continue this simulation.\n");
01539 output_prefix() = output_prefix()+".tmp";
01540 write_restart_file(output_prefix()+".colvars.state");
01541 return cvm::error("Exiting with error until issue is addressed.\n",
01542 COLVARS_INPUT_ERROR);
01543 }
01544
01545 return COLVARS_OK;
01546 }
01547
01548
01549 int colvarmodule::backup_file(char const *filename)
01550 {
01551 return proxy->backup_file(filename);
01552 }
01553
01554
01555 int colvarmodule::write_output_files()
01556 {
01557 int error_code = COLVARS_OK;
01558 cvm::increase_depth();
01559 for (std::vector<colvarbias *>::iterator bi = biases.begin();
01560 bi != biases.end();
01561 bi++) {
01562
01563 if ((*bi)->output_freq == 0 ||
01564 cvm::step_relative() == 0 ||
01565 (cvm::step_absolute() % (*bi)->output_freq) != 0) {
01566 error_code |= (*bi)->write_output_files();
01567 }
01568 error_code |= (*bi)->write_state_to_replicas();
01569 }
01570 cvm::decrease_depth();
01571 return error_code;
01572 }
01573
01574
01575 int colvarmodule::read_traj(char const *traj_filename,
01576 long traj_read_begin,
01577 long traj_read_end)
01578 {
01579 cvm::log("Opening trajectory file \""+
01580 std::string(traj_filename)+"\".\n");
01581
01582
01583
01584 std::ifstream traj_is(traj_filename);
01585
01586 while (true) {
01587 while (true) {
01588
01589 std::string line("");
01590
01591 do {
01592 if (!colvarparse::getline_nocomments(traj_is, line)) {
01593 cvm::log("End of file \""+std::string(traj_filename)+
01594 "\" reached, or corrupted file.\n");
01595 traj_is.close();
01596 return false;
01597 }
01598 } while (line.find_first_not_of(colvarparse::white_space) == std::string::npos);
01599
01600 std::istringstream is(line);
01601
01602 if (!(is >> it)) return false;
01603
01604 if ( (it < traj_read_begin) ) {
01605
01606 if ((it % 1000) == 0)
01607 std::cerr << "Skipping trajectory step " << it
01608 << " \r";
01609
01610 continue;
01611
01612 } else {
01613
01614 if ((it % 1000) == 0)
01615 std::cerr << "Reading from trajectory, step = " << it
01616 << " \r";
01617
01618 if ( (traj_read_end > traj_read_begin) &&
01619 (it > traj_read_end) ) {
01620 std::cerr << "\n";
01621 cvm::error("Reached the end of the trajectory, "
01622 "read_end = "+cvm::to_str(traj_read_end)+"\n",
01623 COLVARS_FILE_ERROR);
01624 return COLVARS_ERROR;
01625 }
01626
01627 for (std::vector<colvar *>::iterator cvi = colvars.begin();
01628 cvi != colvars.end();
01629 cvi++) {
01630 if (!(*cvi)->read_traj(is)) {
01631 cvm::error("Error: in reading colvar \""+(*cvi)->name+
01632 "\" from trajectory file \""+
01633 std::string(traj_filename)+"\".\n",
01634 COLVARS_FILE_ERROR);
01635 return COLVARS_ERROR;
01636 }
01637 }
01638
01639 break;
01640 }
01641 }
01642 }
01643 return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
01644 }
01645
01646
01647 std::ostream & colvarmodule::write_restart(std::ostream &os)
01648 {
01649 os.setf(std::ios::scientific, std::ios::floatfield);
01650 os << "configuration {\n"
01651 << " step " << std::setw(it_width)
01652 << it << "\n"
01653 << " dt " << dt() << "\n"
01654 << " version " << std::string(COLVARS_VERSION) << "\n";
01655 if (proxy->units.size() > 0) {
01656 os << " units " << proxy->units << "\n";
01657 }
01658 os << "}\n\n";
01659
01660 int error_code = COLVARS_OK;
01661
01662 cvm::increase_depth();
01663 for (std::vector<colvar *>::iterator cvi = colvars.begin();
01664 cvi != colvars.end();
01665 cvi++) {
01666 (*cvi)->write_state(os);
01667 }
01668
01669 for (std::vector<colvarbias *>::iterator bi = biases.begin();
01670 bi != biases.end();
01671 bi++) {
01672 (*bi)->write_state(os);
01673 }
01674 cvm::decrease_depth();
01675
01676 if (error_code != COLVARS_OK) {
01677
01678 os.setstate(std::ios::failbit);
01679 }
01680
01681 return os;
01682 }
01683
01684
01685 std::ostream & colvarmodule::write_traj_label(std::ostream &os)
01686 {
01687 os.setf(std::ios::scientific, std::ios::floatfield);
01688
01689 os << "# " << cvm::wrap_string("step", cvm::it_width-2)
01690 << " ";
01691
01692 cvm::increase_depth();
01693 for (std::vector<colvar *>::iterator cvi = colvars.begin();
01694 cvi != colvars.end();
01695 cvi++) {
01696 (*cvi)->write_traj_label(os);
01697 }
01698 for (std::vector<colvarbias *>::iterator bi = biases.begin();
01699 bi != biases.end();
01700 bi++) {
01701 (*bi)->write_traj_label(os);
01702 }
01703 os << "\n";
01704
01705 cvm::decrease_depth();
01706 return os;
01707 }
01708
01709
01710 std::ostream & colvarmodule::write_traj(std::ostream &os)
01711 {
01712 os.setf(std::ios::scientific, std::ios::floatfield);
01713
01714 os << std::setw(cvm::it_width) << it
01715 << " ";
01716
01717 cvm::increase_depth();
01718 for (std::vector<colvar *>::iterator cvi = colvars.begin();
01719 cvi != colvars.end();
01720 cvi++) {
01721 (*cvi)->write_traj(os);
01722 }
01723 for (std::vector<colvarbias *>::iterator bi = biases.begin();
01724 bi != biases.end();
01725 bi++) {
01726 (*bi)->write_traj(os);
01727 }
01728 os << "\n";
01729
01730 cvm::decrease_depth();
01731 return os;
01732 }
01733
01734
01735 void colvarmodule::log(std::string const &message, int min_log_level)
01736 {
01737 if (cvm::log_level() < min_log_level) return;
01738
01739 size_t const d = (cvm::main() != NULL) ? depth() : 0;
01740 if (d > 0) {
01741 proxy->log((std::string(2*d, ' '))+message);
01742 } else {
01743 proxy->log(message);
01744 }
01745 }
01746
01747
01748 void colvarmodule::increase_depth()
01749 {
01750 (depth())++;
01751 }
01752
01753
01754 void colvarmodule::decrease_depth()
01755 {
01756 if (depth() > 0) {
01757 (depth())--;
01758 }
01759 }
01760
01761
01762 size_t & colvarmodule::depth()
01763 {
01764
01765 colvarmodule *cv = cvm::main();
01766 if (proxy->smp_enabled() == COLVARS_OK) {
01767 int const nt = proxy->smp_num_threads();
01768 if (int(cv->depth_v.size()) != nt) {
01769 proxy->smp_lock();
01770
01771 if (cv->depth_v.size() > 0) { cv->depth_s = cv->depth_v[0]; }
01772 cv->depth_v.clear();
01773 cv->depth_v.assign(nt, cv->depth_s);
01774 proxy->smp_unlock();
01775 }
01776 return cv->depth_v[proxy->smp_thread_id()];
01777 }
01778 return cv->depth_s;
01779 }
01780
01781
01782 void colvarmodule::set_error_bits(int code)
01783 {
01784 if (code < 0) {
01785 cvm::log("Error: set_error_bits() received negative error code.\n");
01786 return;
01787 }
01788 proxy->smp_lock();
01789 errorCode |= code | COLVARS_ERROR;
01790 proxy->smp_unlock();
01791 }
01792
01793
01794 bool colvarmodule::get_error_bit(int code)
01795 {
01796 return bool(errorCode & code);
01797 }
01798
01799
01800 void colvarmodule::clear_error()
01801 {
01802 proxy->smp_lock();
01803 errorCode = COLVARS_OK;
01804 proxy->smp_unlock();
01805 proxy->clear_error_msgs();
01806 }
01807
01808
01809 int colvarmodule::error(std::string const &message, int code)
01810 {
01811 set_error_bits(code);
01812 proxy->error(message);
01813 return get_error();
01814 }
01815
01816
01817 int cvm::read_index_file(char const *filename)
01818 {
01819 std::istream &is = proxy->input_stream(filename, "index file");
01820
01821 if (!is) {
01822 return COLVARS_FILE_ERROR;
01823 } else {
01824 index_file_names.push_back(std::string(filename));
01825 }
01826
01827 while (is.good()) {
01828 char open, close;
01829 std::string group_name;
01830 int index_of_group = -1;
01831 if ( (is >> open) && (open == '[') &&
01832 (is >> group_name) &&
01833 (is >> close) && (close == ']') ) {
01834 size_t i = 0;
01835 for ( ; i < index_group_names.size(); i++) {
01836 if (index_group_names[i] == group_name) {
01837
01838 index_of_group = i;
01839 }
01840 }
01841 if (index_of_group < 0) {
01842 index_group_names.push_back(group_name);
01843 index_groups.push_back(NULL);
01844 index_of_group = index_groups.size()-1;
01845 }
01846 } else {
01847 return cvm::error("Error: in parsing index file \""+
01848 std::string(filename)+"\".\n",
01849 COLVARS_INPUT_ERROR);
01850 }
01851
01852 std::vector<int> *old_index_group = index_groups[index_of_group];
01853 std::vector<int> *new_index_group = new std::vector<int>();
01854
01855 int atom_number = 1;
01856 std::streampos pos = is.tellg();
01857 while ( (is >> atom_number) && (atom_number > 0) ) {
01858 new_index_group->push_back(atom_number);
01859 pos = is.tellg();
01860 }
01861
01862 if (old_index_group != NULL) {
01863 bool equal = false;
01864 if (new_index_group->size() == old_index_group->size()) {
01865 if (std::equal(new_index_group->begin(), new_index_group->end(),
01866 old_index_group->begin())) {
01867 equal = true;
01868 }
01869 }
01870 if (! equal) {
01871 new_index_group->clear();
01872 delete new_index_group;
01873 new_index_group = NULL;
01874 return cvm::error("Error: the index group \""+group_name+
01875 "\" was redefined.\n", COLVARS_INPUT_ERROR);
01876 } else {
01877 old_index_group->clear();
01878 delete old_index_group;
01879 old_index_group = NULL;
01880 }
01881 }
01882
01883 index_groups[index_of_group] = new_index_group;
01884
01885 is.clear();
01886 is.seekg(pos, std::ios::beg);
01887 std::string delim;
01888 if ( (is >> delim) && (delim == "[") ) {
01889
01890 is.clear();
01891 is.seekg(pos, std::ios::beg);
01892 } else {
01893 break;
01894 }
01895 }
01896
01897 cvm::log("The following index groups are currently defined:\n");
01898 size_t i = 0;
01899 for ( ; i < index_group_names.size(); i++) {
01900 cvm::log(" "+(index_group_names[i])+" ("+
01901 cvm::to_str((index_groups[i])->size())+" atoms)\n");
01902 }
01903
01904 return proxy->close_input_stream(filename);
01905 }
01906
01907
01908 int colvarmodule::reset_index_groups()
01909 {
01910 size_t i = 0;
01911 for ( ; i < index_groups.size(); i++) {
01912 delete index_groups[i];
01913 index_groups[i] = NULL;
01914 }
01915 index_group_names.clear();
01916 index_groups.clear();
01917 index_file_names.clear();
01918 return COLVARS_OK;
01919 }
01920
01921
01922 int cvm::load_atoms(char const *file_name,
01923 cvm::atom_group &atoms,
01924 std::string const &pdb_field,
01925 double pdb_field_value)
01926 {
01927 return proxy->load_atoms(file_name, atoms, pdb_field, pdb_field_value);
01928 }
01929
01930
01931 int cvm::load_coords(char const *file_name,
01932 std::vector<cvm::rvector> *pos,
01933 cvm::atom_group *atoms,
01934 std::string const &pdb_field,
01935 double pdb_field_value)
01936 {
01937 int error_code = COLVARS_OK;
01938
01939 std::string const ext(strlen(file_name) > 4 ?
01940 (file_name + (strlen(file_name) - 4)) :
01941 file_name);
01942
01943 atoms->create_sorted_ids();
01944
01945 std::vector<cvm::rvector> sorted_pos(atoms->size(), cvm::rvector(0.0));
01946
01947
01948 if (colvarparse::to_lower_cppstr(ext) == std::string(".xyz")) {
01949 if (pdb_field.size() > 0) {
01950 return cvm::error("Error: PDB column may not be specified "
01951 "for XYZ coordinate files.\n", COLVARS_INPUT_ERROR);
01952 }
01953
01954 error_code |= cvm::main()->load_coords_xyz(file_name, &sorted_pos, atoms);
01955 } else {
01956
01957 error_code |= proxy->load_coords(file_name,
01958 sorted_pos, atoms->sorted_ids(),
01959 pdb_field, pdb_field_value);
01960 }
01961
01962 std::vector<int> const &map = atoms->sorted_ids_map();
01963 for (size_t i = 0; i < atoms->size(); i++) {
01964 (*pos)[map[i]] = sorted_pos[i];
01965 }
01966
01967 return error_code;
01968 }
01969
01970
01971 int cvm::load_coords_xyz(char const *filename,
01972 std::vector<rvector> *pos,
01973 cvm::atom_group *atoms)
01974 {
01975 std::istream &xyz_is = proxy->input_stream(filename, "XYZ file");
01976 unsigned int natoms;
01977 char symbol[256];
01978 std::string line;
01979 cvm::real x = 0.0, y = 0.0, z = 0.0;
01980
01981 std::string const error_msg("Error: cannot parse XYZ file \""+
01982 std::string(filename)+"\".\n");
01983
01984 if ( ! (xyz_is >> natoms) ) {
01985 return cvm::error(error_msg, COLVARS_INPUT_ERROR);
01986 }
01987
01988 ++xyz_reader_use_count;
01989 if (xyz_reader_use_count < 2) {
01990 cvm::log("Warning: beginning from 2019-11-26 the XYZ file reader assumes Angstrom units.\n");
01991 }
01992
01993 if (xyz_is.good()) {
01994
01995 cvm::getline(xyz_is, line);
01996 cvm::getline(xyz_is, line);
01997 xyz_is.width(255);
01998 } else {
01999 return cvm::error(error_msg, COLVARS_INPUT_ERROR);
02000 }
02001
02002 std::vector<atom_pos>::iterator pos_i = pos->begin();
02003 size_t xyz_natoms = 0;
02004 if (pos->size() != natoms) {
02005 int next = 0;
02006 std::vector<int>::const_iterator index = atoms->sorted_ids().begin();
02007
02008 for ( ; pos_i != pos->end() ; pos_i++, index++) {
02009 while ( next < *index ) {
02010 cvm::getline(xyz_is, line);
02011 next++;
02012 }
02013 if (xyz_is.good()) {
02014 xyz_is >> symbol;
02015 xyz_is >> x >> y >> z;
02016
02017 (*pos_i)[0] = proxy->angstrom_to_internal(x);
02018 (*pos_i)[1] = proxy->angstrom_to_internal(y);
02019 (*pos_i)[2] = proxy->angstrom_to_internal(z);
02020 xyz_natoms++;
02021 } else {
02022 return cvm::error(error_msg, COLVARS_INPUT_ERROR);
02023 }
02024 }
02025
02026 } else {
02027
02028 for ( ; pos_i != pos->end() ; pos_i++) {
02029 if (xyz_is.good()) {
02030 xyz_is >> symbol;
02031 xyz_is >> x >> y >> z;
02032 (*pos_i)[0] = proxy->angstrom_to_internal(x);
02033 (*pos_i)[1] = proxy->angstrom_to_internal(y);
02034 (*pos_i)[2] = proxy->angstrom_to_internal(z);
02035 xyz_natoms++;
02036 } else {
02037 return cvm::error(error_msg, COLVARS_INPUT_ERROR);
02038 }
02039 }
02040 }
02041
02042 if (xyz_natoms != pos->size()) {
02043 return cvm::error("Error: The number of positions read from file \""+
02044 std::string(filename)+"\" does not match the number of "+
02045 "positions required: "+cvm::to_str(xyz_natoms)+" vs. "+
02046 cvm::to_str(pos->size())+".\n", COLVARS_INPUT_ERROR);
02047 }
02048
02049 return proxy->close_input_stream(filename);
02050 }
02051
02052
02053
02054
02055
02056
02057 cvm::real cvm::dt()
02058 {
02059 return proxy->dt();
02060 }
02061
02062
02063 void cvm::request_total_force()
02064 {
02065 proxy->request_total_force(true);
02066 }
02067
02068
02069 cvm::rvector cvm::position_distance(cvm::atom_pos const &pos1,
02070 cvm::atom_pos const &pos2)
02071 {
02072 return proxy->position_distance(pos1, pos2);
02073 }
02074
02075
02076 cvm::real cvm::rand_gaussian(void)
02077 {
02078 return proxy->rand_gaussian();
02079 }
02080
02081
02082 template<typename T> std::string _to_str(T const &x,
02083 size_t width, size_t prec)
02084 {
02085 std::ostringstream os;
02086 if (width) os.width(width);
02087 if (prec) {
02088 os.setf(std::ios::scientific, std::ios::floatfield);
02089 os.precision(prec);
02090 }
02091 os << x;
02092 return os.str();
02093 }
02094
02095
02096 template<typename T> std::string _to_str_vector(std::vector<T> const &x,
02097 size_t width, size_t prec)
02098 {
02099 if (!x.size()) return std::string("");
02100 std::ostringstream os;
02101 if (prec) {
02102 os.setf(std::ios::scientific, std::ios::floatfield);
02103 }
02104 os << "{ ";
02105 if (width) os.width(width);
02106 if (prec) os.precision(prec);
02107 os << x[0];
02108 for (size_t i = 1; i < x.size(); i++) {
02109 os << ", ";
02110 if (width) os.width(width);
02111 if (prec) os.precision(prec);
02112 os << x[i];
02113 }
02114 os << " }";
02115 return os.str();
02116 }
02117
02118
02119
02120 std::string colvarmodule::to_str(std::string const &x)
02121 {
02122 return std::string("\"")+x+std::string("\"");
02123 }
02124
02125 std::string colvarmodule::to_str(char const *x)
02126 {
02127 return std::string("\"")+std::string(x)+std::string("\"");
02128 }
02129
02130 std::string colvarmodule::to_str(bool x)
02131 {
02132 return (x ? "on" : "off");
02133 }
02134
02135 std::string colvarmodule::to_str(int const &x,
02136 size_t width, size_t prec)
02137 {
02138 return _to_str<int>(x, width, prec);
02139 }
02140
02141 std::string colvarmodule::to_str(size_t const &x,
02142 size_t width, size_t prec)
02143 {
02144 return _to_str<size_t>(x, width, prec);
02145 }
02146
02147 std::string colvarmodule::to_str(long int const &x,
02148 size_t width, size_t prec)
02149 {
02150 return _to_str<long int>(x, width, prec);
02151 }
02152
02153 std::string colvarmodule::to_str(step_number const &x,
02154 size_t width, size_t prec)
02155 {
02156 return _to_str<step_number>(x, width, prec);
02157 }
02158
02159 std::string colvarmodule::to_str(cvm::real const &x,
02160 size_t width, size_t prec)
02161 {
02162 return _to_str<cvm::real>(x, width, prec);
02163 }
02164
02165 std::string colvarmodule::to_str(cvm::rvector const &x,
02166 size_t width, size_t prec)
02167 {
02168 return _to_str<cvm::rvector>(x, width, prec);
02169 }
02170
02171 std::string colvarmodule::to_str(cvm::quaternion const &x,
02172 size_t width, size_t prec)
02173 {
02174 return _to_str<cvm::quaternion>(x, width, prec);
02175 }
02176
02177 std::string colvarmodule::to_str(colvarvalue const &x,
02178 size_t width, size_t prec)
02179 {
02180 return _to_str<colvarvalue>(x, width, prec);
02181 }
02182
02183 std::string colvarmodule::to_str(cvm::vector1d<cvm::real> const &x,
02184 size_t width, size_t prec)
02185 {
02186 return _to_str< cvm::vector1d<cvm::real> >(x, width, prec);
02187 }
02188
02189 std::string colvarmodule::to_str(cvm::matrix2d<cvm::real> const &x,
02190 size_t width, size_t prec)
02191 {
02192 return _to_str< cvm::matrix2d<cvm::real> >(x, width, prec);
02193 }
02194
02195
02196 std::string colvarmodule::to_str(std::vector<int> const &x,
02197 size_t width, size_t prec)
02198 {
02199 return _to_str_vector<int>(x, width, prec);
02200 }
02201
02202 std::string colvarmodule::to_str(std::vector<size_t> const &x,
02203 size_t width, size_t prec)
02204 {
02205 return _to_str_vector<size_t>(x, width, prec);
02206 }
02207
02208 std::string colvarmodule::to_str(std::vector<long int> const &x,
02209 size_t width, size_t prec)
02210 {
02211 return _to_str_vector<long int>(x, width, prec);
02212 }
02213
02214 std::string colvarmodule::to_str(std::vector<cvm::real> const &x,
02215 size_t width, size_t prec)
02216 {
02217 return _to_str_vector<cvm::real>(x, width, prec);
02218 }
02219
02220 std::string colvarmodule::to_str(std::vector<cvm::rvector> const &x,
02221 size_t width, size_t prec)
02222 {
02223 return _to_str_vector<cvm::rvector>(x, width, prec);
02224 }
02225
02226 std::string colvarmodule::to_str(std::vector<cvm::quaternion> const &x,
02227 size_t width, size_t prec)
02228 {
02229 return _to_str_vector<cvm::quaternion>(x, width, prec);
02230 }
02231
02232 std::string colvarmodule::to_str(std::vector<colvarvalue> const &x,
02233 size_t width, size_t prec)
02234 {
02235 return _to_str_vector<colvarvalue>(x, width, prec);
02236 }
02237
02238 std::string colvarmodule::to_str(std::vector<std::string> const &x,
02239 size_t width, size_t prec)
02240 {
02241 return _to_str_vector<std::string>(x, width, prec);
02242 }
02243
02244
02245 std::string cvm::wrap_string(std::string const &s, size_t nchars)
02246 {
02247 if (!s.size()) {
02248 return std::string(nchars, ' ');
02249 } else {
02250 return ( (s.size() <= nchars) ?
02251 (s+std::string(nchars-s.size(), ' ')) :
02252 (std::string(s, 0, nchars)) );
02253 }
02254 }
02255
02256
02257
02258 int colvarmodule::cite_feature(std::string const &feature)
02259 {
02260 return usage_->cite_feature(feature);
02261 }
02262
02263 std::string colvarmodule::feature_report(int flag)
02264 {
02265 return usage_->report(flag);
02266 }
02267
02268
02269 colvarmodule::usage::usage()
02270 {
02271 #include "colvarmodule_refs.h"
02272 }
02273
02274 int colvarmodule::usage::cite_feature(std::string const &feature)
02275 {
02276 if (feature_count_.count(feature) > 0) {
02277 feature_count_[feature] += 1;
02278 return cite_paper(feature_paper_map_[feature]);
02279 }
02280 cvm::log("Warning: cannot cite unknown feature \""+feature+"\"\n");
02281 return COLVARS_OK;
02282 }
02283
02284 int colvarmodule::usage::cite_paper(std::string const &paper)
02285 {
02286 if (paper_count_.count(paper) > 0) {
02287 paper_count_[paper] += 1;
02288 return COLVARS_OK;
02289 }
02290 cvm::log("Warning: cannot cite unknown paper \""+paper+"\"\n");
02291 return COLVARS_OK;
02292 }
02293
02294 std::string colvarmodule::usage::report(int flag)
02295 {
02296 std::string result;
02297 if (flag == 0) {
02298
02299 result += "SUMMARY OF COLVARS FEATURES USED SO FAR AND THEIR CITATIONS:\n";
02300 }
02301 if (flag == 1) {
02302
02303 result += "Colvars module (Fiorin2013, plus other works listed for specific features)\n\n";
02304 }
02305
02306 std::map<std::string, int>::iterator p_iter = paper_count_.begin();
02307 for ( ; p_iter != paper_count_.end(); p_iter++) {
02308 std::string const paper = p_iter->first;
02309 int const count = p_iter->second;
02310 if (count > 0) {
02311 result += "\n";
02312 std::map<std::string, std::string>::iterator f_iter =
02313 feature_paper_map_.begin();
02314 for ( ; f_iter != feature_paper_map_.end(); f_iter++) {
02315 if ((f_iter->second == paper) &&
02316 (feature_count_[f_iter->first] > 0)) {
02317 if (flag == 0) {
02318
02319 result += "- " + f_iter->first + ":\n";
02320 }
02321 if (flag == 1) {
02322
02323 result += "% " + f_iter->first + ":\n";
02324 }
02325 }
02326 }
02327 if (flag == 0) {
02328 result += " " + paper + " " + paper_url_[paper] + "\n";
02329 }
02330 if (flag == 1) {
02331 result += paper_bibtex_[paper] + "\n";
02332 }
02333 }
02334 }
02335
02336 return result;
02337 }
02338
02339
02340
02341 colvarproxy *colvarmodule::proxy = NULL;
02342
02343
02344 cvm::real colvarmodule::debug_gradients_step_size = 1.0e-07;
02345 int colvarmodule::errorCode = 0;
02346 int colvarmodule::log_level_ = 10;
02347 cvm::step_number colvarmodule::it = 0;
02348 cvm::step_number colvarmodule::it_restart = 0;
02349 size_t colvarmodule::restart_out_freq = 0;
02350 size_t colvarmodule::cv_traj_freq = 0;
02351 bool colvarmodule::use_scripted_forces = false;
02352 bool colvarmodule::scripting_after_biases = true;
02353
02354
02355 size_t const colvarmodule::it_width = 12;
02356 size_t const colvarmodule::cv_prec = 14;
02357 size_t const colvarmodule::cv_width = 21;
02358 size_t const colvarmodule::en_prec = 14;
02359 size_t const colvarmodule::en_width = 21;
02360 const char * const colvarmodule::line_marker = (const char *)
02361 "----------------------------------------------------------------------\n";