#include <colvarmodule.h>
Public Types | |
| typedef double | real |
| Defining an abstract real number allows to switch precision. | |
| typedef int | residue_id |
| Residue identifier. | |
| typedef rvector | atom_pos |
| Atom position (different type name from rvector, to make possible future PBC-transparent implementations). | |
| typedef std::vector< atom >::iterator | atom_iter |
| typedef std::vector< atom >::const_iterator | atom_const_iter |
Public Member Functions | |
| colvarmodule (char const *config_name, colvarproxy *proxy_in) | |
| Constructor. | |
| ~colvarmodule () | |
| Destructor. | |
| void | init_colvars (std::string const &conf) |
| Initialize collective variables. | |
| void | init_biases (std::string const &conf) |
| Initialize collective variable biases. | |
| void | change_configuration (std::string const &name, std::string const &conf) |
| Load new configuration - force constant and/or centers only. | |
| real | energy_difference (std::string const &name, std::string const &conf) |
| Calculate change in energy from using alternate configuration. | |
| void | calc () |
| Calculate collective variables and biases. | |
| std::istream & | read_restart (std::istream &is) |
| Read the input restart file. | |
| std::ostream & | write_restart (std::ostream &os) |
| Write the output restart file. | |
| void | write_output_files () |
| Write all output files (called by the proxy). | |
| void | analyze () |
| Perform analysis. | |
| bool | read_traj (char const *traj_filename, size_t traj_read_begin, size_t traj_read_end) |
| Read a collective variable trajectory (post-processing only, not called at runtime). | |
Static Public Member Functions | |
| size_t | step_relative () |
| Return the current step number from the beginning of this run. | |
| size_t | step_absolute () |
| bool | debug () |
| Whether debug output should be enabled (compile-time option). | |
| void | backup_file (char const *filename) |
| Call colvarproxy::backup_file(). | |
| colvar * | colvar_p (std::string const &name) |
| Get the pointer of a colvar from its name (returns NULL if not found). | |
| template<typename T> | |
| std::string | to_str (T const &x, size_t const &width=0, size_t const &prec=0) |
| Quick conversion of an object to a string. | |
| template<typename T> | |
| std::string | to_str (std::vector< T > const &x, size_t const &width=0, size_t const &prec=0) |
| Quick conversion of a vector of objects to a string. | |
| std::string | wrap_string (std::string const &s, size_t const &nchars) |
| Reduce the number of characters in a string. | |
| real | unit_angstrom () |
| Value of the unit for atomic coordinates with respect to angstroms (used by some variables for hard-coded default values). | |
| real | boltzmann () |
| Boltmann constant. | |
| real | temperature () |
| Temperature of the simulation (K). | |
| real | dt () |
| Time step of MD integrator (fs). | |
| void | request_system_force () |
| Request calculation of system force from MD engine. | |
| void | log (std::string const &message) |
| Print a message to the main log. | |
| void | fatal_error (std::string const &message) |
| Print a message to the main log and exit with error code. | |
| void | exit (std::string const &message) |
| Print a message to the main log and exit normally. | |
| rvector | position_distance (atom_pos const &pos1, atom_pos const &pos2) |
| Get the distance between two atomic positions with pbcs handled correctly. | |
| real | position_dist2 (atom_pos const &pos1, atom_pos const &pos2) |
| Get the square distance between two positions (with periodic boundary conditions handled transparently). | |
| void | select_closest_image (atom_pos &pos, atom_pos const &ref_pos) |
| Get the closest periodic image to a reference position. | |
| void | select_closest_images (std::vector< atom_pos > &pos, atom_pos const &ref_pos) |
| Perform select_closest_image() on a set of atomic positions. | |
| void | load_atoms (char const *filename, std::vector< atom > &atoms, std::string const &pdb_field, double const pdb_field_value=0.0) |
| Create atoms from a file. | |
| void | load_coords (char const *filename, std::vector< atom_pos > &pos, const std::vector< int > &indices, std::string const &pdb_field, double const pdb_field_value=0.0) |
| Load the coordinates for a group of atoms from a file (usually a PDB); the number of atoms in "filename" must match the number of elements in "pos". | |
| real | rand_gaussian (void) |
| Pseudo-random number with Gaussian distribution. | |
| void | increase_depth () |
| Increase the depth (number of indentations in the output). | |
| void | decrease_depth () |
| Decrease the depth (number of indentations in the output). | |
Public Attributes | |
| bool | it_restart_from_state_file |
| std::string | restart_out_name |
| Output restart file name. | |
Static Public Attributes | |
| size_t | it = 0 |
| Current step number. | |
| size_t | it_restart = 0 |
| Starting step number for this run. | |
| real | debug_gradients_step_size = 1.0e-03 |
| Finite difference step size (if there is no dynamics, or if gradients need to be tested independently from the size of dt). | |
| std::string | output_prefix = "" |
| Prefix for all output files for this run. | |
| std::string | input_prefix = "" |
| Prefix for files from a previous run (including restart/output). | |
| std::string | restart_in_name = "" |
| input restart file name (determined from input_prefix) | |
| std::vector< colvar * > | colvars |
| Array of collective variables. | |
| std::vector< colvarbias * > | biases |
| Array of collective variable biases. | |
| size_t | n_abf_biases = 0 |
| Number of ABF biases initialized (in normal conditions should be 1). | |
| size_t | n_meta_biases = 0 |
| Number of metadynamics biases initialized (in normal conditions should be 1). | |
| size_t | n_harm_biases = 0 |
| Number of harmonic biases initialized (no limit on the number). | |
| size_t | n_histo_biases = 0 |
| Number of histograms initialized (no limit on the number). | |
| size_t const | it_width = 12 |
| Number of characters to represent a time step. | |
| size_t const | cv_prec = 14 |
| Number of digits to represent a collective variables value(s). | |
| size_t const | cv_width = 21 |
| Number of characters to represent a collective variables value(s). | |
| size_t const | en_prec = 14 |
| Number of digits to represent the collective variables energy. | |
| size_t const | en_width = 21 |
| Number of characters to represent the collective variables energy. | |
| std::string const | line_marker |
| Line separator in the log output. | |
| size_t | cv_traj_freq = 0 |
| Frequency for collective variables trajectory output. | |
| bool | b_analysis = false |
| True if only analysis is performed and not a run. | |
| size_t | restart_out_freq = 0 |
| Frequency for saving output restarts. | |
Protected Attributes | |
| std::ifstream | config_s |
| Configuration file. | |
| colvarparse * | parse |
| Configuration file parser object. | |
| std::string | cv_traj_name |
| Name of the trajectory file. | |
| std::ofstream | cv_traj_os |
| Collective variables output trajectory file. | |
| bool | cv_traj_append |
| Appending to the existing trajectory file? | |
| std::ofstream | restart_out_os |
| Output restart file. | |
Static Protected Attributes | |
| colvarproxy * | proxy = NULL |
| Pointer to the proxy object, used to retrieve atomic data from the hosting program; it is static in order to be accessible from static functions in the colvarmodule class. | |
| size_t | depth = 0 |
| Counter for the current depth in the object hierarchy (useg e.g. in outpu. | |
Friends | |
| class | colvarproxy |
| class | atom |
| class | atom_group |
Class to control the collective variables calculation. An object (usually one) of this class is spawned from the MD program, containing all i/o routines and general interface.
At initialization, the colvarmodule object creates a proxy object to provide a transparent interface between the MD program and the child objects
Definition at line 47 of file colvarmodule.h.
|
|
Definition at line 85 of file colvarmodule.h. |
|
|
Definition at line 84 of file colvarmodule.h. |
|
|
Atom position (different type name from rvector, to make possible future PBC-transparent implementations).
Definition at line 74 of file colvarmodule.h. Referenced by colvar::eigenvector::calc_Jacobian_derivative(), colvar::logmsd::calc_Jacobian_derivative(), colvar::rmsd::calc_Jacobian_derivative(), colvar::dihedral::calc_value(), colvar::angle::calc_value(), colvarmodule::atom_group::center_of_geometry(), colvarmodule::atom_group::center_of_mass(), colvarproxy_namd::load_coords(), colvarmodule::atom_group::parse(), colvarproxy_namd::position_dist2(), position_dist2(), colvarproxy_namd::position_distance(), position_distance(), colvarmodule::atom_group::read_positions(), colvarmodule::atom::reset_data(), colvarproxy_namd::select_closest_image(), select_closest_image(), colvarproxy::select_closest_images(), and select_closest_images(). |
|
|
|
Residue identifier.
Definition at line 61 of file colvarmodule.h. |
|
||||||||||||
|
Constructor.
Definition at line 10 of file colvarmodule.C. References b_analysis, colvarproxy::backup_file(), colvarparse::check_keywords(), COLVARS_VERSION, config_s, cv_traj_append, cv_traj_freq, cv_traj_name, cv_traj_os, debug(), debug_gradients_step_size, fatal_error(), colvarparse::getline_nocomments(), init_biases(), init_colvars(), colvarproxy::input_prefix(), it, it_restart, it_restart_from_state_file, log(), colvarproxy::output_prefix(), output_prefix, parse, proxy, read_restart(), colvarproxy::restart_frequency(), restart_in_name, restart_out_freq, restart_out_name, and colvarproxy::restart_output_prefix(). 00012 {
00013 // pointer to the proxy object
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 // "it_restart" will be set by the input state file, if any;
00027 // "it" should be updated by the proxy
00028 it = it_restart = 0;
00029 it_restart_from_state_file = true;
00030
00031 // open the configfile
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 // read the config file into a string
00038 std::string conf = "";
00039 {
00040 std::string line;
00041 while (colvarparse::getline_nocomments (config_s, line))
00042 conf.append (line+"\n");
00043 // don't need the stream any more
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 // by default overwrite the existing trajectory file
00062 parse->get_keyval (conf, "colvarsTrajAppend", cv_traj_append, false);
00063
00064 // input restart file
00065 restart_in_name = proxy->input_prefix().size() ?
00066 std::string (proxy->input_prefix()+".colvars.state") :
00067 std::string ("") ;
00068
00069 // output restart file
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 // open trajectory file
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 // parse the options for collective variables
00103 init_colvars (conf);
00104
00105 // parse the options for biases
00106 init_biases (conf);
00107
00108 // done with the parsing, check that all keywords are valid
00109 parse->check_keywords (conf, "colvarmodule");
00110 cvm::log (cvm::line_marker);
00111
00112 // read the restart configuration, if available
00113 if (restart_in_name.size()) {
00114 // read the restart file
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 // check if it is possible to save output configuration
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 }
|
|
|
Destructor.
Definition at line 545 of file colvarmodule.C. References biases, colvars, and cv_traj_os. 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 }
|
|
|
Perform analysis.
Definition at line 515 of file colvarmodule.C. References biases, colvars, debug(), decrease_depth(), increase_depth(), it, log(), and step_relative(). 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 // perform colvar-specific analysis
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 // perform bias-specific analysis
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 }
|
|
|
Call colvarproxy::backup_file().
Definition at line 477 of file colvarmodule.h. References colvarproxy::backup_file(), and proxy. Referenced by colvarbias_meta::write_pmf(), and colvarbias_meta::write_replica_state_file(). 00478 {
00479 proxy->backup_file (filename);
00480 }
|
|
|
Boltmann constant.
Definition at line 416 of file colvarmodule.h. References colvarproxy::boltzmann(), and proxy. Referenced by colvar::colvar(), and colvar::update().
|
|
|
Calculate collective variables and biases.
Definition at line 338 of file colvarmodule.C. References colvarproxy::add_energy(), colvarproxy::backup_file(), biases, colvars, cv_traj_freq, cv_traj_name, cv_traj_os, debug(), decrease_depth(), fatal_error(), increase_depth(), it, log(), proxy, restart_out_freq, restart_out_name, restart_out_os, step_absolute(), step_relative(), wrap_string(), and write_restart(). Referenced by colvarproxy_namd::calculate(). 00338 {
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 // calculate collective variables and their gradients
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 // update the biases and communicate their forces to the collective
00360 // variables
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 // sum the forces from all biases for each collective variable
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 // perform runtime analysis of colvars and biases
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 // sum up the forces for each colvar and integrate any internal
00406 // equation of motion
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 // make collective variables communicate their forces to their
00420 // coupled degrees of freedom (i.e. atoms)
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 // write restart file, if needed
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 // write trajectory file, if needed
00448 if (cv_traj_freq) {
00449
00450 if (cvm::debug())
00451 cvm::log ("Writing trajectory file.\n");
00452
00453 // (re)open trajectory file
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 // write labels in the traj file every 1000 lines and at first ts
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 // write collective variable values to the traj file
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 // flush the trajectory file if we are at the restart frequency
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 } // end if (cv_traj_freq)
00512 }
|
|
||||||||||||
|
Load new configuration - force constant and/or centers only.
Definition at line 299 of file colvarmodule.C. References biases, decrease_depth(), fatal_error(), and increase_depth(). 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 }
|
|
|
Get the pointer of a colvar from its name (returns NULL if not found).
Definition at line 513 of file colvar.h. Referenced by colvarbias::add_colvar(), colvar::calc_acf(), and colvar::parse_analysis(). 00514 {
00515 for (std::vector<colvar *>::iterator cvi = cvm::colvars.begin();
00516 cvi != cvm::colvars.end();
00517 cvi++) {
00518 if ((*cvi)->name == name) {
00519 return (*cvi);
00520 }
00521 }
00522 return NULL;
00523 }
|
|
|
|
Decrease the depth (number of indentations in the output).
Definition at line 698 of file colvarmodule.C. References depth. Referenced by analyze(), calc(), colvar::calc(), change_configuration(), colvar::communicate_forces(), energy_difference(), init_biases(), init_colvars(), colvarmodule::atom_group::parse(), read_restart(), colvar::update(), and write_restart(). 00699 {
00700 if (depth) depth--;
00701 }
|
|
|
Time step of MD integrator (fs).
Definition at line 426 of file colvarmodule.h. References colvarproxy::dt(), and proxy. Referenced by colvar::fdiff_velocity(), colvar::update(), and write_restart().
|
|
||||||||||||
|
Calculate change in energy from using alternate configuration.
Definition at line 317 of file colvarmodule.C. References biases, decrease_depth(), fatal_error(), and increase_depth(). 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 }
|
|
|
Print a message to the main log and exit normally.
Definition at line 708 of file colvarmodule.C. References colvarproxy::exit(), and proxy.
|
|
|
|
Increase the depth (number of indentations in the output).
Definition at line 693 of file colvarmodule.C. References depth. Referenced by analyze(), calc(), colvar::calc(), change_configuration(), colvar::communicate_forces(), energy_difference(), init_biases(), init_colvars(), colvarmodule::atom_group::parse(), read_restart(), colvar::update(), and write_restart(). 00694 {
00695 depth++;
00696 }
|
|
|
Initialize collective variable biases. initialize ABF instances initialize harmonic restraints initialize histograms initialize metadynamics instances Definition at line 211 of file colvarmodule.C. References biases, colvarparse::check_keywords(), debug(), decrease_depth(), increase_depth(), colvarparse::key_lookup(), log(), n_abf_biases, n_harm_biases, n_histo_biases, n_meta_biases, and parse. Referenced by colvarmodule(). 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 }
|
|
|
Initialize collective variables.
Definition at line 178 of file colvarmodule.C. References colvarparse::check_keywords(), colvars, debug(), decrease_depth(), fatal_error(), increase_depth(), colvarparse::key_lookup(), log(), and parse. Referenced by colvarmodule(). 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 }
|
|
||||||||||||||||||||
|
Create atoms from a file.
Definition at line 460 of file colvarmodule.h. References colvarproxy::load_atoms(), and proxy. Referenced by colvarmodule::atom_group::parse(). 00464 {
00465 proxy->load_atoms (file_name, atoms, pdb_field, pdb_field_value);
00466 }
|
|
||||||||||||||||||||||||
|
Load the coordinates for a group of atoms from a file (usually a PDB); the number of atoms in "filename" must match the number of elements in "pos".
Definition at line 468 of file colvarmodule.h. References colvarproxy::load_coords(), and proxy. Referenced by colvar::eigenvector::eigenvector(), colvar::orientation::orientation(), and colvarmodule::atom_group::parse(). 00473 {
00474 proxy->load_coords (file_name, pos, indices, pdb_field, pdb_field_value);
00475 }
|
|
|
||||||||||||
|
Get the square distance between two positions (with periodic boundary conditions handled transparently). Note: in the case of periodic boundary conditions, this provides an analytical square distance (while taking the square of position_distance() would produce leads to a cusp) Definition at line 454 of file colvarmodule.h. References atom_pos, colvarproxy::position_dist2(), and proxy. Referenced by colvar::distance_vec::dist2(). 00456 {
00457 return proxy->position_dist2 (pos1, pos2);
00458 }
|
|
||||||||||||
|
Get the distance between two atomic positions with pbcs handled correctly.
Definition at line 448 of file colvarmodule.h. References atom_pos, colvarproxy::position_distance(), and proxy. Referenced by colvar::min_distance::calc_gradients(), colvar::distance_xy::calc_gradients(), colvar::distance_dir::calc_value(), colvar::min_distance::calc_value(), colvar::distance_xy::calc_value(), colvar::distance_z::calc_value(), colvar::distance_vec::calc_value(), colvar::distance::calc_value(), colvar::dihedral::calc_value(), colvar::angle::calc_value(), colvar::distance_vec::dist2_lgrad(), colvar::distance_vec::dist2_rgrad(), colvar::selfcoordnum::switching_function(), and colvar::coordnum::switching_function(). 00450 {
00451 return proxy->position_distance (pos1, pos2);
00452 }
|
|
|
Pseudo-random number with Gaussian distribution.
Definition at line 482 of file colvarmodule.h. References proxy, and colvarproxy::rand_gaussian(). Referenced by colvar::update(). 00483 {
00484 return proxy->rand_gaussian();
00485 }
|
|
|
Read the input restart file.
Definition at line 137 of file colvarmodule.C. References biases, colvars, decrease_depth(), fatal_error(), increase_depth(), it, it_restart, and parse. Referenced by colvarmodule(). 00138 {
00139 {
00140 // read global restart information
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 // colvars restart
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 // biases restart
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 }
|
|
||||||||||||||||
|
Read a collective variable trajectory (post-processing only, not called at runtime).
Definition at line 590 of file colvarmodule.C. References colvars, colvarparse::getline_nocomments(), it, and log(). 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 }
|
|
|
Request calculation of system force from MD engine.
Definition at line 431 of file colvarmodule.h. References proxy, and colvarproxy::request_system_force(). Referenced by colvar::enable(). 00432 {
00433 proxy->request_system_force (true);
00434 }
|
|
||||||||||||
|
Get the closest periodic image to a reference position.
Definition at line 436 of file colvarmodule.h. References atom_pos, proxy, and colvarproxy::select_closest_image(). Referenced by colvarmodule::atom_group::center_of_geometry(), and colvarmodule::atom_group::center_of_mass(). 00438 {
00439 proxy->select_closest_image (pos, ref_pos);
00440 }
|
|
||||||||||||
|
Perform select_closest_image() on a set of atomic positions. After that, distance vectors can then be calculated directly, without using position_distance() Definition at line 442 of file colvarmodule.h. References atom_pos, proxy, and colvarproxy::select_closest_images(). 00444 {
00445 proxy->select_closest_images (pos, ref_pos);
00446 }
|
|
|
Return the current step number from the beginning of the whole calculation Definition at line 100 of file colvarmodule.h. Referenced by calc(), colvarbias_meta::update(), colvarbias_histogram::update(), colvarbias_abf::update(), colvarbias_harmonic::update(), colvarbias_meta::write_replica_state_file(), and colvarbias_meta::write_restart(). 00101 {
00102 return it;
00103 }
|
|
|
Return the current step number from the beginning of this run.
Definition at line 93 of file colvarmodule.h. Referenced by analyze(), calc(), colvar::calc(), colvar::calc_runave(), colvarproxy_namd::calculate(), and colvarbias_abf::update(). 00094 {
00095 return it - it_restart;
00096 }
|
|
|
Temperature of the simulation (K).
Definition at line 421 of file colvarmodule.h. References proxy, and colvarproxy::temperature(). Referenced by colvarbias_abf::colvarbias_abf(), colvar::update(), and colvarbias_meta::write_pmf(). 00422 {
00423 return proxy->temperature();
00424 }
|
|
||||||||||||||||||||
|
Quick conversion of a vector of objects to a string.
Definition at line 385 of file colvarmodule.h. 00387 {
00388 if (!x.size()) return std::string ("");
00389 std::ostringstream os;
00390 if (prec) {
00391 os.setf (std::ios::scientific, std::ios::floatfield);
00392 }
00393 os << "{ ";
00394 if (width) os.width (width);
00395 if (prec) os.precision (prec);
00396 os << x[0];
00397 for (size_t i = 1; i < x.size(); i++) {
00398 os << ", ";
00399 if (width) os.width (width);
00400 if (prec) os.precision (prec);
00401 os << x[i];
00402 }
00403 os << " }";
00404 return os.str();
00405 }
|
|
||||||||||||||||||||
|
Quick conversion of an object to a string.
Definition at line 372 of file colvarmodule.h. Referenced by colvar::alpha_angles::calc_value(), colvarmodule::atom_group::create_sorted_ids(), colvar::dihedPC::dihedPC(), colvarbias_meta::hill::hill(), and colvarmodule::atom_group::parse(). 00374 {
00375 std::ostringstream os;
00376 if (width) os.width (width);
00377 if (prec) {
00378 os.setf (std::ios::scientific, std::ios::floatfield);
00379 os.precision (prec);
00380 }
00381 os << x;
00382 return os.str();
00383 }
|
|
|
Value of the unit for atomic coordinates with respect to angstroms (used by some variables for hard-coded default values).
Definition at line 411 of file colvarmodule.h. References proxy, and colvarproxy::unit_angstrom(). 00412 {
00413 return proxy->unit_angstrom();
00414 }
|
|
||||||||||||
|
Reduce the number of characters in a string.
Definition at line 209 of file colvarmodule.h. Referenced by calc(), colvar::parse_analysis(), and colvar::write_traj_label(). 00211 {
00212 if (!s.size())
00213 return std::string (nchars, ' ');
00214 else
00215 return ( (s.size() <= size_t (nchars)) ?
00216 (s+std::string (nchars-s.size(), ' ')) :
00217 (std::string (s, 0, nchars)) );
00218 }
|
|
|
Write all output files (called by the proxy).
Definition at line 569 of file colvarmodule.C. References colvarproxy::backup_file(), cv_traj_os, log(), output_prefix, proxy, and write_restart(). Referenced by colvarproxy_namd::calculate(). 00570 {
00571 // if this is a simulation run (i.e. not a postprocessing), output data
00572 // must be written to be able to restart the simulation
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 // do not close to avoid problems with multiple NAMD runs
00585 cv_traj_os.flush();
00586 }
|
|
|
Write the output restart file.
Definition at line 658 of file colvarmodule.C. References biases, colvars, decrease_depth(), dt(), increase_depth(), it, and it_width. Referenced by calc(), and write_output_files(). 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 }
|
|
|
Definition at line 82 of file colvarmodule.h. Referenced by colvar::h_bond::h_bond(). |
|
|
Definition at line 83 of file colvarmodule.h. |
|
|
Definition at line 56 of file colvarmodule.h. |
|
|
True if only analysis is performed and not a run.
Definition at line 732 of file colvarmodule.C. Referenced by colvarmodule(). |
|
|
Array of collective variable biases.
Definition at line 717 of file colvarmodule.C. Referenced by analyze(), calc(), change_configuration(), energy_difference(), init_biases(), read_restart(), write_restart(), and ~colvarmodule(). |
|
|
Array of collective variables.
Definition at line 716 of file colvarmodule.C. Referenced by analyze(), calc(), init_colvars(), read_restart(), read_traj(), write_restart(), and ~colvarmodule(). |
|
|
Configuration file.
Definition at line 326 of file colvarmodule.h. Referenced by colvarmodule(). |
|
|
Number of digits to represent a collective variables value(s).
Definition at line 744 of file colvarmodule.C. |
|
|
Appending to the existing trajectory file?
Definition at line 338 of file colvarmodule.h. Referenced by colvarmodule(). |
|
|
Frequency for collective variables trajectory output.
Definition at line 730 of file colvarmodule.C. Referenced by calc(), and colvarmodule(). |
|
|
Name of the trajectory file.
Definition at line 332 of file colvarmodule.h. Referenced by calc(), and colvarmodule(). |
|
|
Collective variables output trajectory file.
Definition at line 335 of file colvarmodule.h. Referenced by calc(), colvarmodule(), write_output_files(), and ~colvarmodule(). |
|
|
Number of characters to represent a collective variables value(s).
Definition at line 745 of file colvarmodule.C. |
|
|
Finite difference step size (if there is no dynamics, or if gradients need to be tested independently from the size of dt).
Definition at line 726 of file colvarmodule.C. Referenced by colvarmodule(). |
|
|
Counter for the current depth in the object hierarchy (useg e.g. in outpu.
Definition at line 731 of file colvarmodule.C. Referenced by decrease_depth(), increase_depth(), and log(). |
|
|
Number of digits to represent the collective variables energy.
Definition at line 746 of file colvarmodule.C. |
|
|
Number of characters to represent the collective variables energy.
Definition at line 747 of file colvarmodule.C. |
|
|
Prefix for files from a previous run (including restart/output).
Definition at line 738 of file colvarmodule.C. |
|
|
Current step number.
Definition at line 727 of file colvarmodule.C. Referenced by analyze(), calc(), colvarproxy_namd::calculate(), colvarmodule(), colvarproxy_namd::colvarproxy_namd(), read_restart(), read_traj(), and write_restart(). |
|
|
Starting step number for this run.
Definition at line 728 of file colvarmodule.C. Referenced by colvarmodule(), colvarproxy_namd::colvarproxy_namd(), and read_restart(). |
|
|
If true, get it_restart from the state file; if set to false, the MD program is providing it Definition at line 107 of file colvarmodule.h. Referenced by colvarmodule(). |
|
|
Number of characters to represent a time step.
Definition at line 743 of file colvarmodule.C. Referenced by write_restart(). |
|
|
Initial value:
"----------------------------------------------------------------------\n"
Definition at line 748 of file colvarmodule.C. |
|
|
Number of ABF biases initialized (in normal conditions should be 1).
Definition at line 718 of file colvarmodule.C. Referenced by init_biases(). |
|
|
Number of harmonic biases initialized (no limit on the number).
Definition at line 719 of file colvarmodule.C. Referenced by init_biases(). |
|
|
Number of histograms initialized (no limit on the number).
Definition at line 720 of file colvarmodule.C. Referenced by init_biases(). |
|
|
Number of metadynamics biases initialized (in normal conditions should be 1).
Definition at line 721 of file colvarmodule.C. Referenced by init_biases(). |
|
|
Prefix for all output files for this run.
Definition at line 737 of file colvarmodule.C. Referenced by colvarmodule(), and write_output_files(). |
|
|
Configuration file parser object.
Definition at line 329 of file colvarmodule.h. Referenced by colvarmodule(), init_biases(), init_colvars(), and read_restart(). |
|
|
Pointer to the proxy object, used to retrieve atomic data from the hosting program; it is static in order to be accessible from static functions in the colvarmodule class.
Definition at line 722 of file colvarmodule.C. Referenced by backup_file(), boltzmann(), calc(), colvarmodule(), dt(), exit(), fatal_error(), load_atoms(), load_coords(), log(), position_dist2(), position_distance(), rand_gaussian(), request_system_force(), select_closest_image(), select_closest_images(), temperature(), unit_angstrom(), and write_output_files(). |
|
|
input restart file name (determined from input_prefix)
Definition at line 739 of file colvarmodule.C. Referenced by colvarmodule(). |
|
|
Frequency for saving output restarts.
Definition at line 729 of file colvarmodule.C. Referenced by calc(), and colvarmodule(). |
|
|
Output restart file name.
Definition at line 319 of file colvarmodule.h. Referenced by calc(), and colvarmodule(). |
|
|
Output restart file.
Definition at line 341 of file colvarmodule.h. Referenced by calc(). |
1.3.9.1