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

colvarmodule.h

Go to the documentation of this file.
00001 #ifndef COLVARMODULE_H
00002 #define COLVARMODULE_H
00003 
00004 #ifndef COLVARS_VERSION
00005 #define COLVARS_VERSION "2013-01-30"
00006 #endif
00007 
00008 #ifndef COLVARS_DEBUG
00009 #define COLVARS_DEBUG false
00010 #endif
00011 
00020 
00021 
00022 #include <iostream>
00023 #include <iomanip>
00024 #include <string>
00025 #include <sstream>
00026 #include <fstream>
00027 #include <cmath>
00028 #include <vector>
00029 #include <list>
00030 
00031 
00032 class colvarparse;
00033 class colvar;
00034 class colvarbias;
00035 class colvarproxy;
00036 
00037 
00047 class colvarmodule {
00048 
00049 private:
00050 
00052   colvarmodule();
00053 
00054 public:
00055 
00056   friend class colvarproxy;
00057   
00059   typedef  double    real;
00061   typedef  int       residue_id;
00062 
00063   class rvector;
00064   template <class T,
00065             size_t const length> class vector1d;
00066   template <class T,
00067             size_t const outer_length,
00068             size_t const inner_length> class matrix2d;
00069   class quaternion;
00070   class rotation;
00071 
00074   typedef rvector atom_pos;
00075 
00077   class rmatrix;
00078 
00079   // allow these classes to access protected data
00080   class atom;
00081   class atom_group;
00082   friend class atom;
00083   friend class atom_group;
00084   typedef std::vector<atom>::iterator       atom_iter;
00085   typedef std::vector<atom>::const_iterator atom_const_iter;
00086 
00088   static size_t it;
00090   static size_t it_restart;
00091 
00093   static inline size_t step_relative()
00094   {
00095     return it - it_restart;
00096   }
00097 
00100   static inline size_t step_absolute()
00101   {
00102     return it;
00103   }
00104 
00107   bool it_restart_from_state_file;
00108 
00112   static real debug_gradients_step_size;
00113 
00115   static std::string output_prefix;
00116 
00118   static std::string input_prefix;
00119 
00121   static std::string restart_in_name;
00122 
00123 
00125   static std::vector<colvar *>     colvars;
00126 
00127   /* TODO: implement named CVCs
00129   static std::vector<cvc *>     cvcs;
00131   inline void register_cvc (cvc *p) {
00132     cvcs.push_back(p);
00133   }
00134   */
00135 
00137   static std::vector<colvarbias *> biases;
00140   static size_t n_abf_biases;
00143   static size_t n_meta_biases;
00146   static size_t n_harm_biases;
00149   static size_t n_histo_biases;
00150    
00152   static inline bool debug()
00153   {
00154     return COLVARS_DEBUG;
00155   }
00156 
00157 
00160   colvarmodule (char const *config_name, 
00161                 colvarproxy *proxy_in);
00163   ~colvarmodule();
00164 
00166   void init_colvars (std::string const &conf);
00167 
00169   void init_biases (std::string const &conf);
00170 
00173   void change_configuration (std::string const &bias_name, std::string const &conf);
00174 
00177   real energy_difference (std::string const &bias_name, std::string const &conf);
00178 
00180   void calc();
00182   std::istream & read_restart (std::istream &is);
00184   std::ostream & write_restart (std::ostream &os);
00186   void write_output_files();
00188   static void backup_file (char const *filename);
00189 
00191   void analyze();
00194   bool read_traj (char const *traj_filename,
00195                   size_t      traj_read_begin,
00196                   size_t      traj_read_end);
00197  
00199   static colvar * colvar_p (std::string const &name);
00200 
00202   template<typename T> static std::string to_str (T const &x, 
00203                                                   size_t const &width = 0,
00204                                                   size_t const &prec = 0);
00206   template<typename T> static std::string to_str (std::vector<T> const &x, 
00207                                                   size_t const &width = 0,
00208                                                   size_t const &prec = 0);
00209 
00211   static inline std::string wrap_string (std::string const &s,
00212                                          size_t const &nchars)
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   }
00221 
00223   static size_t const it_width;
00225   static size_t const cv_prec;
00227   static size_t const cv_width;
00229   static size_t const en_prec;
00231   static size_t const en_width;
00233   static std::string const line_marker;
00234 
00235 
00236   // proxy functions
00237 
00240   static real unit_angstrom();
00241 
00243   static real boltzmann();
00244 
00246   static real temperature();
00247 
00249   static real dt();
00250   
00252   static void request_system_force();
00253   
00255   static void log (std::string const &message);
00256 
00258   static void fatal_error (std::string const &message);
00259 
00261   static void exit (std::string const &message);
00262 
00263 
00266   static rvector position_distance (atom_pos const &pos1,
00267                                     atom_pos const &pos2);
00268 
00269 
00276   static real position_dist2 (atom_pos const &pos1,
00277                               atom_pos const &pos2);
00278 
00282   static void select_closest_image (atom_pos &pos,
00283                                     atom_pos const &ref_pos);
00284 
00289   static void select_closest_images (std::vector<atom_pos> &pos,
00290                                      atom_pos const &ref_pos);
00291 
00292 
00297   static void load_atoms (char const *filename,
00298                           std::vector<atom> &atoms,
00299                           std::string const &pdb_field,
00300                           double const pdb_field_value = 0.0);
00301 
00305   static void load_coords (char const *filename,
00306                            std::vector<atom_pos> &pos,
00307                            const std::vector<int> &indices,
00308                            std::string const &pdb_field,
00309                            double const pdb_field_value = 0.0);
00310 
00311 
00313   static size_t cv_traj_freq;
00314 
00316   static bool   b_analysis;
00317 
00319   static size_t restart_out_freq;
00321   std::string   restart_out_name;
00322 
00324   static real rand_gaussian (void);
00325 protected:
00326 
00328   std::ifstream config_s;
00329 
00331   colvarparse *parse;
00332 
00334   std::string   cv_traj_name;
00335 
00337   std::ofstream cv_traj_os;
00338 
00340   bool          cv_traj_append;
00341 
00343   std::ofstream restart_out_os;
00344 
00348   static colvarproxy *proxy;
00349 
00351   static size_t depth;
00352 
00353 public:
00354 
00356   static void increase_depth();
00357 
00359   static void decrease_depth();
00360 };
00361 
00362 
00364 typedef colvarmodule cvm;
00365 
00366 
00367 #include "colvartypes.h"
00368 
00369 
00370 std::ostream & operator << (std::ostream &os, cvm::rvector const &v);
00371 std::istream & operator >> (std::istream &is, cvm::rvector &v);
00372 
00373 
00374 template<typename T> std::string cvm::to_str (T const &x, 
00375                                               size_t const &width,
00376                                               size_t const &prec) {
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 }
00386 
00387 template<typename T> std::string cvm::to_str (std::vector<T> const &x, 
00388                                               size_t const &width,
00389                                               size_t const &prec) {
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 }
00408 
00409 
00410 #include "colvarproxy.h"
00411 
00412 
00413 inline cvm::real cvm::unit_angstrom()
00414 {
00415   return proxy->unit_angstrom();
00416 }
00417 
00418 inline cvm::real cvm::boltzmann()
00419 {
00420   return proxy->boltzmann();
00421 }
00422 
00423 inline cvm::real cvm::temperature()
00424 {
00425   return proxy->temperature();
00426 }
00427 
00428 inline cvm::real cvm::dt()
00429 {
00430   return proxy->dt();
00431 }
00432 
00433 inline void cvm::request_system_force()
00434 {
00435   proxy->request_system_force (true);
00436 }
00437   
00438 inline void cvm::select_closest_image (atom_pos &pos,
00439                                        atom_pos const &ref_pos)
00440 {
00441   proxy->select_closest_image (pos, ref_pos);
00442 }
00443 
00444 inline void cvm::select_closest_images (std::vector<atom_pos> &pos,
00445                                         atom_pos const &ref_pos)
00446 {
00447   proxy->select_closest_images (pos, ref_pos);
00448 }
00449 
00450 inline cvm::rvector cvm::position_distance (atom_pos const &pos1,
00451                                             atom_pos const &pos2)
00452 {
00453   return proxy->position_distance (pos1, pos2);
00454 }
00455 
00456 inline cvm::real cvm::position_dist2 (cvm::atom_pos const &pos1,
00457                                       cvm::atom_pos const &pos2)
00458 {
00459   return proxy->position_dist2 (pos1, pos2);
00460 }
00461 
00462 inline void cvm::load_atoms (char const *file_name,
00463                              std::vector<cvm::atom> &atoms,
00464                              std::string const &pdb_field,
00465                              double const pdb_field_value)
00466 {
00467   proxy->load_atoms (file_name, atoms, pdb_field, pdb_field_value);
00468 }
00469 
00470 inline void cvm::load_coords (char const *file_name,
00471                               std::vector<cvm::atom_pos> &pos,
00472                               const std::vector<int> &indices,
00473                               std::string const &pdb_field,
00474                               double const pdb_field_value)
00475 {
00476   proxy->load_coords (file_name, pos, indices, pdb_field, pdb_field_value);
00477 }
00478 
00479 inline void cvm::backup_file (char const *filename)
00480 {
00481   proxy->backup_file (filename);
00482 }
00483 
00484 inline cvm::real cvm::rand_gaussian (void)
00485 {
00486   return proxy->rand_gaussian();
00487 }
00488 
00489 #endif
00490 
00491 
00492 // Emacs
00493 // Local Variables:
00494 // mode: C++
00495 // End:

Generated on Sun May 19 04:07:44 2013 for NAMD by  doxygen 1.3.9.1