Main Page | Namespace List | Class Hierarchy | Alphabetical List | Class List | File List | Namespace Members | Class Members | File Members

colvarmodule Class Reference

Collective variables module (main class). More...

#include <colvarmodule.h>

List of all members.

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 &bias_name, std::string const &conf)
real energy_difference (std::string const &bias_name, std::string const &conf)
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().
colvarcolvar_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.
colvarparseparse
 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

colvarproxyproxy = 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


Detailed Description

Collective variables module (main class).

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.


Member Typedef Documentation

typedef std::vector<atom>::const_iterator colvarmodule::atom_const_iter
 

Definition at line 85 of file colvarmodule.h.

typedef std::vector<atom>::iterator colvarmodule::atom_iter
 

Definition at line 84 of file colvarmodule.h.

typedef rvector colvarmodule::atom_pos
 

Atom position (different type name from rvector, to make possible future PBC-transparent implementations).

Definition at line 74 of file colvarmodule.h.

Referenced by colvarmodule::atom_group::calc_apply_roto_translation(), colvarmodule::atom_group::calc_fit_gradients(), colvar::eigenvector::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(), colvarmodule::atom_group::center_ref_pos(), colvarproxy_namd::load_coords(), colvarmodule::atom_group::parse(), colvarproxy_namd::position_dist2(), position_dist2(), colvarproxy_namd::position_distance(), position_distance(), colvarmodule::atom::reset_data(), colvarproxy_namd::select_closest_image(), select_closest_image(), colvarproxy::select_closest_images(), and select_closest_images().

typedef double colvarmodule::real
 

Defining an abstract real number allows to switch precision.

Definition at line 59 of file colvarmodule.h.

Referenced by colvar::alpha_angles::apply_force(), colvar_grid_gradient::average(), colvar::calc(), colvarmodule::atom_group::calc_fit_gradients(), colvar::rmsd::calc_gradients(), colvar::gyration::calc_gradients(), colvar::distance_inv::calc_gradients(), colvar::gyration::calc_Jacobian_derivative(), colvar::calc_runave(), colvar::alpha_angles::calc_value(), colvar::rmsd::calc_value(), colvar::distance_inv::calc_value(), colvar::cvc::compare(), colvar::orientation::orientation(), colvarproxy_namd::position_dist2(), colvarbias_meta::project_hills(), colvarvalue::reset(), colvarbias_abf::update(), colvarbias_harmonic::update(), colvar::update(), colvar_grid_gradient::value_output(), colvar_grid_scalar::value_output(), colvar_grid_gradient::write_1D_integral(), and colvar::write_acf().

typedef int colvarmodule::residue_id
 

Residue identifier.

Definition at line 61 of file colvarmodule.h.


Constructor & Destructor Documentation

colvarmodule::colvarmodule char const *  config_name,
colvarproxy proxy_in
 

Constructor.

Parameters:
config_name Configuration file name
restart_name (optional) Restart file name

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_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   parse->get_keyval (conf, "debugGradientsStepSize", debug_gradients_step_size, 1.0e-07,
00050                      colvarparse::parse_silent);
00051 
00052   parse->get_keyval (conf, "eigenvalueCrossingThreshold",
00053                      colvarmodule::rotation::crossing_threshold, 1.0e-02,
00054                      colvarparse::parse_silent);
00055 
00056   parse->get_keyval (conf, "colvarsTrajFrequency", cv_traj_freq, 100);
00057   parse->get_keyval (conf, "colvarsRestartFrequency", restart_out_freq,
00058                      proxy->restart_frequency());
00059 
00060   // by default overwrite the existing trajectory file
00061   parse->get_keyval (conf, "colvarsTrajAppend", cv_traj_append, false);
00062 
00063   // input restart file
00064   restart_in_name = proxy->input_prefix().size() ?
00065     std::string (proxy->input_prefix()+".colvars.state") :
00066     std::string ("") ;
00067 
00068   // output restart file
00069   restart_out_name = proxy->restart_output_prefix().size() ?
00070     std::string (proxy->restart_output_prefix()+".colvars.state") :
00071     std::string ("");
00072 
00073   if (restart_out_name.size())
00074     cvm::log ("The restart output state file will be \""+restart_out_name+"\".\n");
00075 
00076   output_prefix = proxy->output_prefix();
00077 
00078   cvm::log ("The final output state file will be \""+
00079             (output_prefix.size() ?
00080              std::string (output_prefix+".colvars.state") :
00081              std::string ("colvars.state"))+"\".\n");
00082 
00083   cv_traj_name = 
00084     (output_prefix.size() ?
00085      std::string (output_prefix+".colvars.traj") :
00086      std::string ("colvars.traj"));
00087   cvm::log ("The trajectory file will be \""+
00088             cv_traj_name+"\".\n");
00089 
00090   if (cv_traj_freq) {
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 
00103   // parse the options for collective variables
00104   init_colvars (conf);
00105 
00106   // parse the options for biases
00107   init_biases (conf);
00108 
00109   // done with the parsing, check that all keywords are valid
00110   parse->check_keywords (conf, "colvarmodule");
00111   cvm::log (cvm::line_marker);
00112 
00113   // read the restart configuration, if available
00114   if (restart_in_name.size()) {
00115     // read the restart file
00116     std::ifstream input_is (restart_in_name.c_str());
00117     if (!input_is.good())
00118       fatal_error ("Error: in opening restart file \""+
00119                    std::string (restart_in_name)+"\".\n");
00120     else {
00121       cvm::log ("Restarting from file \""+restart_in_name+"\".\n");
00122       read_restart (input_is);
00123       cvm::log (cvm::line_marker);
00124     }
00125   }
00126 
00127   // check if it is possible to save output configuration
00128   if ((!output_prefix.size()) && (!restart_out_name.size())) {
00129     cvm::fatal_error ("Error: neither the final output state file or "
00130                       "the output restart file could be defined, exiting.\n");
00131   }
00132 
00133   cvm::log ("Collective variables module initialized.\n");
00134   cvm::log (cvm::line_marker);
00135 }

colvarmodule::~colvarmodule  ) 
 

Destructor.

Definition at line 546 of file colvarmodule.C.

References biases, colvars, and cv_traj_os.

00547 {
00548   for (std::vector<colvar *>::iterator cvi = colvars.begin();
00549        cvi != colvars.end();
00550        cvi++) {
00551     delete *cvi;
00552   }
00553   colvars.clear();
00554 
00555   for (std::vector<colvarbias *>::iterator bi = biases.begin();
00556        bi != biases.end();
00557        bi++) {
00558     delete *bi;
00559   }
00560   biases.clear();
00561 
00562   if (cv_traj_os.good()) {
00563     cv_traj_os.close();
00564   }
00565 
00566   delete parse;
00567 }  


Member Function Documentation

void colvarmodule::analyze  ) 
 

Perform analysis.

Definition at line 516 of file colvarmodule.C.

References biases, colvars, debug(), decrease_depth(), increase_depth(), it, log(), and step_relative().

00517 {
00518   if (cvm::debug()) {
00519     cvm::log ("colvarmodule::analyze(), step = "+cvm::to_str (it)+".\n");
00520   }
00521 
00522   if (cvm::step_relative() == 0)
00523     cvm::log ("Performing analysis.\n");
00524 
00525   // perform colvar-specific analysis
00526   for (std::vector<colvar *>::iterator cvi = colvars.begin();
00527        cvi != colvars.end();
00528        cvi++) {
00529     cvm::increase_depth();
00530     (*cvi)->analyse(); 
00531     cvm::decrease_depth();
00532   }
00533 
00534   // perform bias-specific analysis
00535   for (std::vector<colvarbias *>::iterator bi = biases.begin();
00536        bi != biases.end();
00537        bi++) {
00538     cvm::increase_depth();
00539     (*bi)->analyse(); 
00540     cvm::decrease_depth();
00541   }
00542 
00543 }

void cvm::backup_file char const *  filename  )  [inline, static]
 

Call colvarproxy::backup_file().

Definition at line 479 of file colvarmodule.h.

References colvarproxy::backup_file(), and proxy.

Referenced by colvarbias_meta::write_pmf(), and colvarbias_meta::write_replica_state_file().

00480 {
00481   proxy->backup_file (filename);
00482 }

cvm::real cvm::boltzmann  )  [inline, static]
 

Boltmann constant.

Definition at line 418 of file colvarmodule.h.

References colvarproxy::boltzmann(), and proxy.

Referenced by colvar::colvar(), and colvar::update().

00419 {
00420   return proxy->boltzmann();
00421 }

void colvarmodule::calc  ) 
 

Calculate collective variables and biases.

Definition at line 339 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().

00339                         {
00340   cvm::real total_bias_energy = 0.0;
00341   cvm::real total_colvar_energy = 0.0;
00342 
00343   if (cvm::debug()) {
00344     cvm::log (cvm::line_marker);
00345     cvm::log ("Collective variables module, step no. "+
00346               cvm::to_str (cvm::step_absolute())+"\n");
00347   }
00348 
00349   // calculate collective variables and their gradients
00350   if (cvm::debug()) 
00351     cvm::log ("Calculating collective variables.\n");
00352   cvm::increase_depth();
00353   for (std::vector<colvar *>::iterator cvi = colvars.begin();
00354        cvi != colvars.end();
00355        cvi++) {
00356     (*cvi)->calc(); 
00357   }
00358   cvm::decrease_depth();
00359 
00360   // update the biases and communicate their forces to the collective
00361   // variables
00362   if (cvm::debug() && biases.size())
00363     cvm::log ("Updating collective variable biases.\n");
00364   cvm::increase_depth();
00365   for (std::vector<colvarbias *>::iterator bi = biases.begin();
00366        bi != biases.end();
00367        bi++) {
00368     total_bias_energy += (*bi)->update(); 
00369   }
00370   cvm::decrease_depth();
00371 
00372   // sum the forces from all biases for each collective variable
00373   if (cvm::debug() && biases.size())
00374     cvm::log ("Collecting forces from all biases.\n");
00375   cvm::increase_depth();
00376   for (std::vector<colvar *>::iterator cvi = colvars.begin();
00377        cvi != colvars.end();
00378        cvi++) {
00379     (*cvi)->reset_bias_force(); 
00380   }
00381   for (std::vector<colvarbias *>::iterator bi = biases.begin();
00382        bi != biases.end();
00383        bi++) {
00384     (*bi)->communicate_forces(); 
00385   }
00386   cvm::decrease_depth();
00387 
00388   if (cvm::b_analysis) {
00389     // perform runtime analysis of colvars and biases
00390     if (cvm::debug() && biases.size())
00391       cvm::log ("Perform runtime analyses.\n");
00392     cvm::increase_depth();
00393     for (std::vector<colvar *>::iterator cvi = colvars.begin();
00394          cvi != colvars.end();
00395          cvi++) {
00396       (*cvi)->analyse(); 
00397     }
00398     for (std::vector<colvarbias *>::iterator bi = biases.begin();
00399          bi != biases.end();
00400          bi++) {
00401       (*bi)->analyse(); 
00402     }
00403     cvm::decrease_depth();
00404   }
00405 
00406   // sum up the forces for each colvar and integrate any internal
00407   // equation of motion
00408   if (cvm::debug())
00409     cvm::log ("Updating the internal degrees of freedom "
00410               "of colvars (if they have any).\n");
00411   cvm::increase_depth();
00412   for (std::vector<colvar *>::iterator cvi = colvars.begin();
00413        cvi != colvars.end();
00414        cvi++) {
00415     total_colvar_energy += (*cvi)->update();
00416   }
00417   cvm::decrease_depth();
00418   proxy->add_energy (total_bias_energy + total_colvar_energy);
00419 
00420   // make collective variables communicate their forces to their
00421   // coupled degrees of freedom (i.e. atoms)
00422   if (cvm::debug())
00423     cvm::log ("Communicating forces from the colvars to the atoms.\n");
00424   cvm::increase_depth();
00425   for (std::vector<colvar *>::iterator cvi = colvars.begin();
00426        cvi != colvars.end();
00427        cvi++) {
00428     if ((*cvi)->tasks[colvar::task_gradients])
00429       (*cvi)->communicate_forces(); 
00430   }
00431   cvm::decrease_depth();
00432 
00433   // write restart file, if needed
00434   if (restart_out_freq && restart_out_name.size() && !cvm::b_analysis) {
00435     if ( (cvm::step_relative() > 0) &&
00436          ((cvm::step_absolute() % restart_out_freq) == 0) ) {
00437       cvm::log ("Writing the state file \""+
00438                 restart_out_name+"\".\n");
00439       proxy->backup_file (restart_out_name.c_str());
00440       restart_out_os.open (restart_out_name.c_str());
00441       restart_out_os.setf (std::ios::scientific, std::ios::floatfield);
00442       if (!write_restart (restart_out_os))
00443         cvm::fatal_error ("Error: in writing restart file.\n");
00444       restart_out_os.close();
00445     }    
00446   }
00447 
00448   // write trajectory file, if needed
00449   if (cv_traj_freq) {
00450 
00451     if (cvm::debug())
00452       cvm::log ("Writing trajectory file.\n");
00453 
00454     // (re)open trajectory file
00455     if (!cv_traj_os.good()) {
00456       if (cv_traj_append) {
00457         cvm::log ("Appending to colvar trajectory file \""+cv_traj_name+
00458                   "\".\n");
00459         cv_traj_os.open (cv_traj_name.c_str(), std::ios::app);
00460       } else {
00461         cvm::log ("Overwriting colvar trajectory file \""+cv_traj_name+
00462                   "\".\n");
00463         proxy->backup_file (cv_traj_name.c_str());
00464         cv_traj_os.open (cv_traj_name.c_str(), std::ios::out);
00465       }
00466       cv_traj_os.setf (std::ios::scientific, std::ios::floatfield);
00467     }
00468 
00469     // write labels in the traj file every 1000 lines and at first ts
00470     cvm::increase_depth();
00471     if ((cvm::step_absolute() % (cv_traj_freq * 1000)) == 0 || cvm::step_relative() == 0) {
00472       cv_traj_os << "# " << cvm::wrap_string ("step", cvm::it_width-2)
00473                  << " ";
00474       if (cvm::debug()) 
00475         cv_traj_os.flush();
00476       for (std::vector<colvar *>::iterator cvi = colvars.begin();
00477            cvi != colvars.end();
00478            cvi++) {
00479         (*cvi)->write_traj_label (cv_traj_os);
00480       }
00481       cv_traj_os << "\n";
00482       if (cvm::debug()) 
00483         cv_traj_os.flush();
00484     }
00485     cvm::decrease_depth();
00486 
00487     // write collective variable values to the traj file
00488     cvm::increase_depth();
00489     if ((cvm::step_absolute() % cv_traj_freq) == 0) {
00490       cv_traj_os << std::setw (cvm::it_width) << it
00491                  << " ";
00492       for (std::vector<colvar *>::iterator cvi = colvars.begin();
00493            cvi != colvars.end();
00494            cvi++) {
00495         (*cvi)->write_traj (cv_traj_os);
00496       }
00497       cv_traj_os << "\n";
00498       if (cvm::debug())
00499         cv_traj_os.flush();
00500     }
00501     cvm::decrease_depth();
00502 
00503     if (restart_out_freq) { 
00504       // flush the trajectory file if we are at the restart frequency
00505       if ( (cvm::step_relative() > 0) &&
00506            ((cvm::step_absolute() % restart_out_freq) == 0) ) {
00507         cvm::log ("Synchronizing (emptying the buffer of) trajectory file \""+
00508                   cv_traj_name+"\".\n");
00509         cv_traj_os.flush();
00510       }
00511     }
00512   } // end if (cv_traj_freq)
00513 }

void colvarmodule::change_configuration std::string const &  bias_name,
std::string const &  conf
 

Load new configuration for the given bias - currently works for harmonic (force constant and/or centers)

Definition at line 300 of file colvarmodule.C.

References biases, decrease_depth(), fatal_error(), and increase_depth().

00302 {
00303   cvm::increase_depth();
00304   int found = 0;
00305   for (std::vector<colvarbias *>::iterator bi = biases.begin();
00306        bi != biases.end();
00307        bi++) {
00308     if ( (*bi)->name == bias_name ) {
00309       ++found;
00310       (*bi)->change_configuration (conf);
00311     }
00312   }
00313   if (found < 1) cvm::fatal_error ("Error: bias not found.\n");
00314   if (found > 1) cvm::fatal_error ("Error: duplicate bias name.\n");
00315   cvm::decrease_depth();
00316 }

colvar * cvm::colvar_p std::string const &  name  )  [inline, static]
 

Get the pointer of a colvar from its name (returns NULL if not found).

Definition at line 524 of file colvar.h.

Referenced by colvarbias::add_colvar(), colvar::calc_acf(), and colvar::parse_analysis().

00525 {
00526   for (std::vector<colvar *>::iterator cvi = cvm::colvars.begin();
00527        cvi != cvm::colvars.end();
00528        cvi++) {
00529     if ((*cvi)->name == name) {
00530       return (*cvi);
00531     }
00532   }
00533   return NULL;
00534 }

bool colvarmodule::debug  )  [inline, static]
 

Whether debug output should be enabled (compile-time option).

Definition at line 152 of file colvarmodule.h.

Referenced by colvarbias::add_colvar(), colvar_grid< size_t >::address(), colvar::alpha_angles::alpha_angles(), analyze(), colvarmodule::atom::atom(), colvar::build_atom_list(), calc(), colvar::calc(), colvarmodule::atom_group::calc_fit_gradients(), colvar::calc_runave(), colvar::alpha_angles::calc_value(), colvarproxy_namd::calculate(), colvarparse::check_keywords(), colvar::colvar(), colvar_grid< size_t >::colvar_grid(), colvarbias_harmonic::colvarbias_harmonic(), colvarbias_meta::colvarbias_meta(), colvarproxy_namd::colvarproxy_namd(), colvarbias::communicate_forces(), colvar::communicate_forces(), colvarbias_meta::create_hill(), colvar::cvc::cvc(), colvarbias_meta::delete_hill(), colvar::dihedPC::dihedPC(), colvar::dihedral::dihedral(), colvar::enable(), colvarproxy_namd::fatal_error(), colvar::h_bond::h_bond(), colvarbias_meta::hill::hill(), init_biases(), init_colvars(), colvar_grid< size_t >::map_grid(), colvarmodule::atom_group::parse(), colvarbias_meta::project_hills(), colvarbias_meta::read_hill(), colvarbias_meta::read_replica_files(), colvarbias_meta::read_restart(), colvarbias_meta::update(), colvarbias_histogram::update(), colvarbias_abf::update(), colvarbias_harmonic::update(), colvar::update(), and colvarbias_meta::update_replicas_registry().

00153   {
00154     return COLVARS_DEBUG;
00155   }

void cvm::decrease_depth  )  [static]
 

Decrease the depth (number of indentations in the output).

Definition at line 707 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(), write_output_files(), and write_restart().

00708 {
00709   if (depth) depth--;
00710 }

cvm::real cvm::dt  )  [inline, static]
 

Time step of MD integrator (fs).

Definition at line 428 of file colvarmodule.h.

References colvarproxy::dt(), and proxy.

Referenced by colvar::fdiff_velocity(), colvar::update(), and write_restart().

00429 {
00430   return proxy->dt();
00431 }

cvm::real colvarmodule::energy_difference std::string const &  bias_name,
std::string const &  conf
 

Calculate change in energy from using alt. config. for the given bias - currently works for harmonic (force constant and/or centers)

Definition at line 318 of file colvarmodule.C.

References biases, decrease_depth(), fatal_error(), and increase_depth().

00320 {
00321   cvm::increase_depth();
00322   cvm::real energy_diff = 0.;
00323   int found = 0;
00324   for (std::vector<colvarbias *>::iterator bi = biases.begin();
00325        bi != biases.end();
00326        bi++) {
00327     if ( (*bi)->name == bias_name ) {
00328       ++found;
00329       energy_diff = (*bi)->energy_difference (conf);
00330     }
00331   }
00332   if (found < 1) cvm::fatal_error ("Error: bias not found.\n");
00333   if (found > 1) cvm::fatal_error ("Error: duplicate bias name.\n");
00334   cvm::decrease_depth();
00335   return energy_diff;
00336 }

void cvm::exit std::string const &  message  )  [static]
 

Print a message to the main log and exit normally.

Definition at line 717 of file colvarmodule.C.

References colvarproxy::exit(), and proxy.

00718 {
00719   proxy->exit (message);
00720 }

void cvm::fatal_error std::string const &  message  )  [static]
 

Print a message to the main log and exit with error code.

Definition at line 712 of file colvarmodule.C.

References colvarproxy::fatal_error(), and proxy.

Referenced by colvarmodule::atom_group::add_atom(), colvarbias::add_colvar(), colvar_grid< size_t >::add_grid(), colvar_grid< size_t >::address(), colvar::alpha_angles::alpha_angles(), colvarmodule::atom_group::apply_colvar_force(), colvarmodule::atom_group::apply_force(), colvarmodule::atom_group::apply_forces(), colvarmodule::atom::atom(), calc(), colvar::calc_acf(), colvar::cvc::calc_force_invgrads(), colvar::cvc::calc_Jacobian_derivative(), colvarproxy_namd::calculate(), change_configuration(), colvarbias::change_configuration(), colvar_grid< size_t >::check_consistency(), colvarparse::check_keywords(), colvarvalue::check_types(), colvar::colvar(), colvar_grid< size_t >::colvar_grid(), colvarbias::colvarbias(), colvarbias_abf::colvarbias_abf(), colvarbias_harmonic::colvarbias_harmonic(), colvarbias_histogram::colvarbias_histogram(), colvarbias_meta::colvarbias_meta(), colvarmodule(), colvar::orientation::compare(), colvar::distance_dir::compare(), colvar::distance_vec::compare(), colvar::cvc::compare(), colvar::coordnum::coordnum(), colvar_grid< size_t >::create(), colvarmodule::atom_group::create_sorted_ids(), colvar::dihedPC::dihedPC(), colvar::distance_inv::distance_inv(), colvar::distance_z::distance_z(), colvar::eigenvector::eigenvector(), colvar::enable(), energy_difference(), colvarbias::energy_difference(), colvarvalue::error_lside(), colvarvalue::error_rside(), colvar_grid_scalar::gradient_finite_diff(), colvar::h_bond::h_bond(), colvar::inertia_z::inertia_z(), init_colvars(), jacobi(), colvarparse::key_lookup(), colvarproxy_namd::load_atoms(), colvarproxy_namd::load_coords(), colvar_grid< size_t >::map_grid(), operator>>(), colvarmodule::quaternion::operator[](), colvar::orientation::orientation(), colvarvalue::p2leg_opt(), colvarmodule::atom_group::parse(), colvar::parse_analysis(), colvar::cvc::parse_group(), colvar_grid< size_t >::parse_params(), pdb_field_str2enum(), colvar::periodic_boundaries(), colvarmodule::atom_group::positions(), colvarmodule::atom_group::positions_shifted(), colvarbias_meta::read_hill(), colvar_grid< size_t >::read_multicol(), colvar_grid< size_t >::read_raw_error(), read_restart(), colvarbias_meta::read_restart(), colvarbias_harmonic::read_restart(), colvar::read_restart(), colvarmodule::atom::read_velocity(), colvar::rmsd::rmsd(), colvar::selfcoordnum::selfcoordnum(), colvarmodule::atom_group::system_force(), colvarmodule::atom_group::system_forces(), colvarvalue::undef_op(), colvarbias_meta::update(), colvarbias_histogram::update(), colvarbias_meta::update_replicas_registry(), colvar_grid_scalar::value_input(), colvar_grid_scalar::value_output(), colvarmodule::atom_group::velocities(), colvar_grid< size_t >::wrap(), colvar_grid_gradient::write_1D_integral(), colvar::write_output_files(), and colvarbias_meta::write_replica_state_file().

00713 {
00714   proxy->fatal_error (message);
00715 }

void cvm::increase_depth  )  [static]
 

Increase the depth (number of indentations in the output).

Definition at line 702 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(), write_output_files(), and write_restart().

00703 {
00704   depth++;
00705 }

void colvarmodule::init_biases std::string const &  conf  ) 
 

Initialize collective variable biases.

initialize ABF instances

initialize harmonic restraints

initialize histograms

initialize metadynamics instances

Definition at line 212 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().

00213 {
00214   if (cvm::debug())
00215     cvm::log ("Initializing the collective variables biases.\n");
00216 
00217   {
00219     std::string abf_conf = "";
00220     size_t abf_pos = 0;
00221     while (parse->key_lookup (conf, "abf", abf_conf, abf_pos)) {
00222       if (abf_conf.size()) {
00223         cvm::log (cvm::line_marker);
00224         cvm::increase_depth();
00225         biases.push_back (new colvarbias_abf (abf_conf, "abf"));
00226         (biases.back())->check_keywords (abf_conf, "abf");
00227         cvm::decrease_depth();
00228         n_abf_biases++;
00229       } else {
00230         cvm::log ("Warning: \"abf\" keyword found without configuration.\n");
00231       }
00232       abf_conf = "";
00233     }
00234   }
00235 
00236   {
00238     std::string harm_conf = "";
00239     size_t harm_pos = 0;
00240     while (parse->key_lookup (conf, "harmonic", harm_conf, harm_pos)) {
00241       if (harm_conf.size()) {
00242         cvm::log (cvm::line_marker);
00243         cvm::increase_depth();
00244         biases.push_back (new colvarbias_harmonic (harm_conf, "harmonic"));
00245         (biases.back())->check_keywords (harm_conf, "harmonic");
00246         cvm::decrease_depth();
00247         n_harm_biases++;
00248       } else {
00249         cvm::log ("Warning: \"harmonic\" keyword found without configuration.\n");
00250       }
00251       harm_conf = "";
00252     }
00253   }
00254 
00255   {
00257     std::string histo_conf = "";
00258     size_t histo_pos = 0;
00259     while (parse->key_lookup (conf, "histogram", histo_conf, histo_pos)) {
00260       if (histo_conf.size()) {
00261         cvm::log (cvm::line_marker);
00262         cvm::increase_depth();
00263         biases.push_back (new colvarbias_histogram (histo_conf, "histogram"));
00264         (biases.back())->check_keywords (histo_conf, "histogram");
00265         cvm::decrease_depth();
00266         n_histo_biases++;
00267       } else {
00268         cvm::log ("Warning: \"histogram\" keyword found without configuration.\n");
00269       }
00270       histo_conf = "";
00271     }
00272   }
00273 
00274   {
00276     std::string meta_conf = "";
00277     size_t meta_pos = 0;
00278     while (parse->key_lookup (conf, "metadynamics", meta_conf, meta_pos)) {
00279       if (meta_conf.size()) {
00280         cvm::log (cvm::line_marker);
00281         cvm::increase_depth();
00282         biases.push_back (new colvarbias_meta (meta_conf, "metadynamics"));
00283         (biases.back())->check_keywords (meta_conf, "metadynamics");
00284         cvm::decrease_depth();
00285         n_meta_biases++;
00286       } else {
00287         cvm::log ("Warning: \"metadynamics\" keyword found without configuration.\n");
00288       }
00289       meta_conf = "";
00290     }
00291   }
00292 
00293   if (biases.size())
00294     cvm::log (cvm::line_marker);
00295   cvm::log ("Collective variables biases initialized, "+
00296             cvm::to_str (biases.size())+" in total.\n");
00297 }

void colvarmodule::init_colvars std::string const &  conf  ) 
 

Initialize collective variables.

Definition at line 179 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().

00180 {
00181   if (cvm::debug())
00182     cvm::log ("Initializing the collective variables.\n");
00183 
00184   std::string colvar_conf = "";
00185   size_t pos = 0;
00186   while (parse->key_lookup (conf, "colvar", colvar_conf, pos)) {
00187 
00188     if (colvar_conf.size()) {
00189       cvm::log (cvm::line_marker);
00190       cvm::increase_depth();
00191       colvars.push_back (new colvar (colvar_conf));
00192       (colvars.back())->check_keywords (colvar_conf, "colvar");
00193       cvm::decrease_depth();
00194     } else {
00195       cvm::log ("Warning: \"colvar\" keyword found without any configuration.\n"); 
00196     }
00197     colvar_conf = "";
00198   }
00199 
00200 
00201   if (!colvars.size())
00202     cvm::fatal_error ("Error: no collective variables defined.\n");
00203 
00204   if (colvars.size())
00205     cvm::log (cvm::line_marker);
00206   cvm::log ("Collective variables initialized, "+
00207             cvm::to_str (colvars.size())+
00208             " in total.\n");
00209 }

void cvm::load_atoms char const *  filename,
std::vector< atom > &  atoms,
std::string const &  pdb_field,
double const   pdb_field_value = 0.0
[inline, static]
 

Create atoms from a file.

Parameters:
filename name of the file (usually a PDB)
atoms array of the atoms to be allocated
pdb_field (optiona) if "filename" is a PDB file, use this field to determine which are the atoms to be set

Definition at line 462 of file colvarmodule.h.

References colvarproxy::load_atoms(), and proxy.

Referenced by colvarmodule::atom_group::parse().

00466 {
00467   proxy->load_atoms (file_name, atoms, pdb_field, pdb_field_value);
00468 }

void cvm::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
[inline, static]
 

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 470 of file colvarmodule.h.

References colvarproxy::load_coords(), and proxy.

Referenced by colvar::eigenvector::eigenvector(), colvar::orientation::orientation(), colvarmodule::atom_group::parse(), and colvar::rmsd::rmsd().

00475 {
00476   proxy->load_coords (file_name, pos, indices, pdb_field, pdb_field_value);
00477 }

void cvm::log std::string const &  message  )  [static]
 

Print a message to the main log.

Definition at line 694 of file colvarmodule.C.

References depth, colvarproxy::log(), and proxy.

Referenced by colvarbias::add_colvar(), colvar::alpha_angles::alpha_angles(), analyze(), colvar::angle::angle(), colvarmodule::atom::atom(), colvarmodule::atom_group::atom_group(), colvar::build_atom_list(), calc(), colvar::calc(), colvar::calc_acf(), colvarmodule::atom_group::calc_fit_gradients(), colvar::tilt::calc_gradients(), colvar::eigenvector::calc_gradients(), colvar::rmsd::calc_gradients(), colvar::inertia_z::calc_gradients(), colvar::inertia::calc_gradients(), colvar::gyration::calc_gradients(), colvar::calc_runave(), colvar::alpha_angles::calc_value(), colvarproxy_namd::calculate(), colvarparse::check_keywords(), colvarvalue::check_types(), colvar::colvar(), colvar_grid< size_t >::colvar_grid(), colvarbias::colvarbias(), colvarbias_abf::colvarbias_abf(), colvarbias_harmonic::colvarbias_harmonic(), colvarbias_histogram::colvarbias_histogram(), colvarbias_meta::colvarbias_meta(), colvarmodule(), colvarproxy_namd::colvarproxy_namd(), colvarbias::communicate_forces(), colvar::communicate_forces(), colvar::cvc::cvc(), colvar::cvc::debug_gradients(), colvarbias_meta::delete_hill(), colvar::dihedPC::dihedPC(), colvar::dihedral::dihedral(), colvar::distance::distance(), colvar::distance_z::distance_z(), colvar::eigenvector::eigenvector(), colvar::enable(), colvarproxy_namd::exit(), colvarproxy_namd::fatal_error(), colvar::gyration::gyration(), colvar::h_bond::h_bond(), colvarbias_meta::hill::hill(), colvar::inertia_z::inertia_z(), init_biases(), init_colvars(), colvar_grid< size_t >::map_grid(), colvar::orientation::orientation(), colvarmodule::atom_group::parse(), colvar::parse_analysis(), colvarbias_meta::project_hills(), colvarbias_meta::read_hill(), colvar_grid< size_t >::read_multicol(), colvarbias_meta::read_replica_files(), colvar_grid_gradient::read_restart(), colvar_grid_scalar::read_restart(), colvar_grid_count::read_restart(), colvarbias_meta::read_restart(), colvarbias_harmonic::read_restart(), colvar::read_restart(), read_traj(), colvar::read_traj(), colvar::rmsd::rmsd(), colvar::spin_angle::spin_angle(), colvar::tilt::tilt(), colvarbias_meta::update(), colvarbias_histogram::update(), colvarbias_abf::update(), colvarbias_harmonic::update(), colvar::update(), colvarbias_meta::update_replicas_registry(), colvar::write_acf(), write_output_files(), and colvar::write_output_files().

00695 {
00696   if (depth > 0)
00697     proxy->log ((std::string (2*depth, ' '))+message);
00698   else
00699     proxy->log (message);
00700 }

cvm::real cvm::position_dist2 atom_pos const &  pos1,
atom_pos const &  pos2
[inline, static]
 

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 456 of file colvarmodule.h.

References atom_pos, colvarproxy::position_dist2(), and proxy.

Referenced by colvar::distance_vec::dist2().

00458 {
00459   return proxy->position_dist2 (pos1, pos2);
00460 }

cvm::rvector cvm::position_distance atom_pos const &  pos1,
atom_pos const &  pos2
[inline, static]
 

Get the distance between two atomic positions with pbcs handled correctly.

Definition at line 450 of file colvarmodule.h.

References atom_pos, colvarproxy::position_distance(), and proxy.

Referenced by colvar::distance_xy::calc_gradients(), colvar::distance_inv::calc_value(), colvar::distance_dir::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().

00452 {
00453   return proxy->position_distance (pos1, pos2);
00454 }

cvm::real cvm::rand_gaussian void   )  [inline, static]
 

Pseudo-random number with Gaussian distribution.

Definition at line 484 of file colvarmodule.h.

References proxy, and colvarproxy::rand_gaussian().

Referenced by colvar::update().

00485 {
00486   return proxy->rand_gaussian();
00487 }

std::istream & colvarmodule::read_restart std::istream &  is  ) 
 

Read the input restart file.

Definition at line 138 of file colvarmodule.C.

References biases, colvars, decrease_depth(), fatal_error(), increase_depth(), it, it_restart, and parse.

Referenced by colvarmodule().

00139 {
00140   {
00141     // read global restart information
00142     std::string restart_conf;
00143     if (is >> colvarparse::read_block ("configuration", restart_conf)) {
00144       if (it_restart_from_state_file) {
00145         parse->get_keyval (restart_conf, "step",
00146                            it_restart, (size_t) 0,
00147                            colvarparse::parse_silent);
00148         it = it_restart;
00149       }
00150     }
00151     is.clear();
00152   }
00153 
00154   // colvars restart
00155   cvm::increase_depth();
00156   for (std::vector<colvar *>::iterator cvi = colvars.begin();
00157        cvi != colvars.end();
00158        cvi++) {
00159     if ( !((*cvi)->read_restart (is)) )
00160       cvm::fatal_error ("Error: in reading restart configuration for collective variable \""+
00161                         (*cvi)->name+"\".\n");
00162   }
00163 
00164   // biases restart
00165   for (std::vector<colvarbias *>::iterator bi = biases.begin();
00166        bi != biases.end();
00167        bi++) {
00168     if (!((*bi)->read_restart (is)))
00169       fatal_error ("Error: in reading restart configuration for bias \""+
00170                    (*bi)->name+"\".\n");
00171   }
00172   cvm::decrease_depth();
00173 
00174   return is;
00175 }

bool colvarmodule::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).

Definition at line 599 of file colvarmodule.C.

References colvars, colvarparse::getline_nocomments(), it, and log().

00602 {
00603   cvm::log ("Opening trajectory file \""+
00604             std::string (traj_filename)+"\".\n");
00605   std::ifstream traj_is (traj_filename);
00606 
00607   while (true) {
00608     while (true) {
00609 
00610       std::string line ("");
00611 
00612       do {
00613         if (!colvarparse::getline_nocomments (traj_is, line)) {
00614           cvm::log ("End of file \""+std::string (traj_filename)+
00615                     "\" reached, or corrupted file.\n");
00616           traj_is.close();
00617           return false;
00618         }
00619       } while (line.find_first_not_of (colvarparse::white_space) == std::string::npos);
00620 
00621       std::istringstream is (line);
00622 
00623       if (!(is >> it)) return false;
00624 
00625       if ( (it < traj_read_begin) ) {
00626 
00627         if ((it % 1000) == 0)
00628           std::cerr << "Skipping trajectory step " << it
00629                     << "                    \r";
00630 
00631         continue;
00632 
00633       } else { 
00634 
00635         if ((it % 1000) == 0)
00636           std::cerr << "Reading from trajectory, step = " << it
00637                     << "                    \r";
00638 
00639         if ( (traj_read_end > traj_read_begin) &&
00640              (it > traj_read_end) ) {
00641           std::cerr << "\n";
00642           cvm::log ("Reached the end of the trajectory, "
00643                     "read_end = "+cvm::to_str (traj_read_end)+"\n");
00644           return false;
00645         }
00646 
00647         for (std::vector<colvar *>::iterator cvi = colvars.begin();
00648              cvi != colvars.end();
00649              cvi++) {
00650           if (!(*cvi)->read_traj (is)) {
00651             cvm::log ("Error: in reading colvar \""+(*cvi)->name+
00652                       "\" from trajectory file \""+
00653                       std::string (traj_filename)+"\".\n");
00654             return false;
00655           }
00656         }
00657 
00658         break;
00659       }
00660     }
00661   }
00662 
00663   return true;
00664 }

void cvm::request_system_force  )  [inline, static]
 

Request calculation of system force from MD engine.

Definition at line 433 of file colvarmodule.h.

References proxy, and colvarproxy::request_system_force().

Referenced by colvar::enable().

00434 {
00435   proxy->request_system_force (true);
00436 }

void cvm::select_closest_image atom_pos pos,
atom_pos const &  ref_pos
[inline, static]
 

Get the closest periodic image to a reference position.

Parameters:
pos The position to look for the closest periodic image
ref_pos (optional) The reference position

Definition at line 438 of file colvarmodule.h.

References atom_pos, proxy, and colvarproxy::select_closest_image().

00440 {
00441   proxy->select_closest_image (pos, ref_pos);
00442 }

void cvm::select_closest_images std::vector< atom_pos > &  pos,
atom_pos const &  ref_pos
[inline, static]
 

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 444 of file colvarmodule.h.

References atom_pos, proxy, and colvarproxy::select_closest_images().

00446 {
00447   proxy->select_closest_images (pos, ref_pos);
00448 }

size_t colvarmodule::step_absolute  )  [inline, static]
 

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   }

size_t colvarmodule::step_relative  )  [inline, static]
 

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(), colvarbias_abf::update(), and colvarbias_harmonic::update().

00094   {
00095     return it - it_restart;
00096   }

cvm::real cvm::temperature  )  [inline, static]
 

Temperature of the simulation (K).

Definition at line 423 of file colvarmodule.h.

References proxy, and colvarproxy::temperature().

Referenced by colvarbias_abf::colvarbias_abf(), colvar::update(), and colvarbias_meta::write_pmf().

00424 {
00425   return proxy->temperature();
00426 }

template<typename T>
std::string cvm::to_str std::vector< T > const &  x,
size_t const &  width = 0,
size_t const &  prec = 0
[static]
 

Quick conversion of a vector of objects to a string.

Definition at line 387 of file colvarmodule.h.

00389                                                                   {
00390   if (!x.size()) return std::string ("");
00391   std::ostringstream os;
00392   if (prec) {
00393     os.setf (std::ios::scientific, std::ios::floatfield);
00394   }
00395   os << "{ ";
00396   if (width) os.width (width);
00397   if (prec) os.precision (prec);
00398   os << x[0];
00399   for (size_t i = 1; i < x.size(); i++) {
00400     os << ", ";
00401     if (width) os.width (width);
00402     if (prec) os.precision (prec);
00403     os << x[i];
00404   }
00405   os << " }";
00406   return os.str();
00407 }

template<typename T>
std::string cvm::to_str T const &  x,
size_t const &  width = 0,
size_t const &  prec = 0
[static]
 

Quick conversion of an object to a string.

Definition at line 374 of file colvarmodule.h.

Referenced by colvar::alpha_angles::calc_value(), colvarmodule::atom_group::create_sorted_ids(), colvar::dihedPC::dihedPC(), colvarbias_meta::hill::hill(), colvarmodule::atom_group::parse(), and colvar::rmsd::rmsd().

00376                                                                   {
00377   std::ostringstream os;
00378   if (width) os.width (width);
00379   if (prec) {
00380     os.setf (std::ios::scientific, std::ios::floatfield);
00381     os.precision (prec);
00382   }
00383   os << x;
00384   return os.str();
00385 }

cvm::real cvm::unit_angstrom  )  [inline, static]
 

Value of the unit for atomic coordinates with respect to angstroms (used by some variables for hard-coded default values).

Definition at line 413 of file colvarmodule.h.

References proxy, and colvarproxy::unit_angstrom().

00414 {
00415   return proxy->unit_angstrom();
00416 }

std::string colvarmodule::wrap_string std::string const &  s,
size_t const &  nchars
[inline, static]
 

Reduce the number of characters in a string.

Definition at line 211 of file colvarmodule.h.

Referenced by calc(), colvar::parse_analysis(), and colvar::write_traj_label().

00213   {
00214     if (!s.size())
00215       return std::string (nchars, ' ');
00216     else
00217       return ( (s.size() <= size_t (nchars)) ?
00218                (s+std::string (nchars-s.size(), ' ')) :
00219                (std::string (s, 0, nchars)) );
00220   }

void colvarmodule::write_output_files  ) 
 

Write all output files (called by the proxy).

Definition at line 570 of file colvarmodule.C.

References colvarproxy::backup_file(), colvars, cv_traj_os, decrease_depth(), increase_depth(), log(), output_prefix, proxy, and write_restart().

Referenced by colvarproxy_namd::calculate().

00571 {
00572   // if this is a simulation run (i.e. not a postprocessing), output data
00573   // must be written to be able to restart the simulation
00574   std::string const out_name =
00575     (output_prefix.size() ?
00576      std::string (output_prefix+".colvars.state") :
00577      std::string ("colvars.state"));
00578   cvm::log ("Saving collective variables state to \""+out_name+"\".\n");
00579   proxy->backup_file (out_name.c_str());
00580   std::ofstream out (out_name.c_str());
00581   out.setf (std::ios::scientific, std::ios::floatfield);
00582   this->write_restart (out);
00583   out.close();
00584 
00585   cvm::increase_depth();
00586   for (std::vector<colvar *>::iterator cvi = colvars.begin();
00587        cvi != colvars.end();
00588        cvi++) {
00589     (*cvi)->write_output_files();
00590   }
00591   cvm::decrease_depth();
00592   
00593   // do not close to avoid problems with multiple NAMD runs
00594   cv_traj_os.flush();
00595 }

std::ostream & colvarmodule::write_restart std::ostream &  os  ) 
 

Write the output restart file.

Definition at line 667 of file colvarmodule.C.

References biases, colvars, decrease_depth(), dt(), increase_depth(), it, and it_width.

Referenced by calc(), and write_output_files().

00668 {
00669   os << "configuration {\n"
00670      << "  step " << std::setw (it_width)
00671      << it << "\n"
00672      << "  dt " << dt() << "\n"
00673      << "}\n\n";
00674 
00675   cvm::increase_depth();
00676   for (std::vector<colvar *>::iterator cvi = colvars.begin();
00677        cvi != colvars.end();
00678        cvi++) {
00679     (*cvi)->write_restart (os);
00680   }
00681 
00682   for (std::vector<colvarbias *>::iterator bi = biases.begin();
00683        bi != biases.end();
00684        bi++) {
00685     (*bi)->write_restart (os);
00686   }
00687   cvm::decrease_depth();
00688 
00689   return os;
00690 }


Friends And Related Function Documentation

friend class atom [friend]
 

Definition at line 82 of file colvarmodule.h.

Referenced by colvar::h_bond::h_bond().

friend class atom_group [friend]
 

Definition at line 83 of file colvarmodule.h.

friend class colvarproxy [friend]
 

Definition at line 56 of file colvarmodule.h.


Member Data Documentation

bool colvarmodule::b_analysis = false [static]
 

True if only analysis is performed and not a run.

Definition at line 741 of file colvarmodule.C.

Referenced by colvarmodule().

std::vector< colvarbias * > colvarmodule::biases [static]
 

Array of collective variable biases.

Definition at line 726 of file colvarmodule.C.

Referenced by analyze(), calc(), change_configuration(), energy_difference(), init_biases(), read_restart(), write_restart(), and ~colvarmodule().

std::vector< colvar * > colvarmodule::colvars [static]
 

Array of collective variables.

Definition at line 725 of file colvarmodule.C.

Referenced by analyze(), calc(), init_colvars(), read_restart(), read_traj(), write_output_files(), write_restart(), and ~colvarmodule().

std::ifstream colvarmodule::config_s [protected]
 

Configuration file.

Definition at line 328 of file colvarmodule.h.

Referenced by colvarmodule().

size_t const colvarmodule::cv_prec = 14 [static]
 

Number of digits to represent a collective variables value(s).

Definition at line 753 of file colvarmodule.C.

bool colvarmodule::cv_traj_append [protected]
 

Appending to the existing trajectory file?

Definition at line 340 of file colvarmodule.h.

Referenced by colvarmodule().

size_t colvarmodule::cv_traj_freq = 0 [static]
 

Frequency for collective variables trajectory output.

Definition at line 739 of file colvarmodule.C.

Referenced by calc(), and colvarmodule().

std::string colvarmodule::cv_traj_name [protected]
 

Name of the trajectory file.

Definition at line 334 of file colvarmodule.h.

Referenced by calc(), and colvarmodule().

std::ofstream colvarmodule::cv_traj_os [protected]
 

Collective variables output trajectory file.

Definition at line 337 of file colvarmodule.h.

Referenced by calc(), colvarmodule(), write_output_files(), and ~colvarmodule().

size_t const colvarmodule::cv_width = 21 [static]
 

Number of characters to represent a collective variables value(s).

Definition at line 754 of file colvarmodule.C.

cvm::real colvarmodule::debug_gradients_step_size = 1.0e-03 [static]
 

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 735 of file colvarmodule.C.

Referenced by colvarmodule().

size_t colvarmodule::depth = 0 [static, protected]
 

Counter for the current depth in the object hierarchy (useg e.g. in outpu.

Definition at line 740 of file colvarmodule.C.

Referenced by decrease_depth(), increase_depth(), and log().

size_t const colvarmodule::en_prec = 14 [static]
 

Number of digits to represent the collective variables energy.

Definition at line 755 of file colvarmodule.C.

size_t const colvarmodule::en_width = 21 [static]
 

Number of characters to represent the collective variables energy.

Definition at line 756 of file colvarmodule.C.

std::string colvarmodule::input_prefix = "" [static]
 

Prefix for files from a previous run (including restart/output).

Definition at line 747 of file colvarmodule.C.

size_t colvarmodule::it = 0 [static]
 

Current step number.

Definition at line 736 of file colvarmodule.C.

Referenced by analyze(), calc(), colvarproxy_namd::calculate(), colvarmodule(), colvarproxy_namd::colvarproxy_namd(), read_restart(), read_traj(), and write_restart().

size_t colvarmodule::it_restart = 0 [static]
 

Starting step number for this run.

Definition at line 737 of file colvarmodule.C.

Referenced by colvarmodule(), colvarproxy_namd::colvarproxy_namd(), and read_restart().

bool colvarmodule::it_restart_from_state_file
 

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().

size_t const colvarmodule::it_width = 12 [static]
 

Number of characters to represent a time step.

Definition at line 752 of file colvarmodule.C.

Referenced by write_restart().

std::string const colvarmodule::line_marker [static]
 

Initial value:

  "----------------------------------------------------------------------\n"
Line separator in the log output.

Definition at line 757 of file colvarmodule.C.

size_t colvarmodule::n_abf_biases = 0 [static]
 

Number of ABF biases initialized (in normal conditions should be 1).

Definition at line 727 of file colvarmodule.C.

Referenced by init_biases().

size_t colvarmodule::n_harm_biases = 0 [static]
 

Number of harmonic biases initialized (no limit on the number).

Definition at line 728 of file colvarmodule.C.

Referenced by init_biases().

size_t colvarmodule::n_histo_biases = 0 [static]
 

Number of histograms initialized (no limit on the number).

Definition at line 729 of file colvarmodule.C.

Referenced by init_biases().

size_t colvarmodule::n_meta_biases = 0 [static]
 

Number of metadynamics biases initialized (in normal conditions should be 1).

Definition at line 730 of file colvarmodule.C.

Referenced by init_biases().

std::string colvarmodule::output_prefix = "" [static]
 

Prefix for all output files for this run.

Definition at line 746 of file colvarmodule.C.

Referenced by colvarmodule(), and write_output_files().

colvarparse* colvarmodule::parse [protected]
 

Configuration file parser object.

Definition at line 331 of file colvarmodule.h.

Referenced by colvarmodule(), init_biases(), init_colvars(), and read_restart().

colvarproxy * colvarmodule::proxy = NULL [static, protected]
 

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 731 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().

std::string colvarmodule::restart_in_name = "" [static]
 

input restart file name (determined from input_prefix)

Definition at line 748 of file colvarmodule.C.

Referenced by colvarmodule().

size_t colvarmodule::restart_out_freq = 0 [static]
 

Frequency for saving output restarts.

Definition at line 738 of file colvarmodule.C.

Referenced by calc(), and colvarmodule().

std::string colvarmodule::restart_out_name
 

Output restart file name.

Definition at line 321 of file colvarmodule.h.

Referenced by calc(), and colvarmodule().

std::ofstream colvarmodule::restart_out_os [protected]
 

Output restart file.

Definition at line 343 of file colvarmodule.h.

Referenced by calc().


The documentation for this class was generated from the following files:
Generated on Wed May 22 04:07:22 2013 for NAMD by  doxygen 1.3.9.1