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

colvar.h

Go to the documentation of this file.
00001 // -*- c++ -*-
00002 
00003 #ifndef COLVAR_H
00004 #define COLVAR_H
00005 
00006 #include <iostream>
00007 #include <iomanip>
00008 #include <cmath>
00009 
00010 #include "colvarmodule.h"
00011 #include "colvarvalue.h"
00012 #include "colvarparse.h"
00013 
00014 
00040 
00041 class colvar : public colvarparse {
00042 
00043 public:
00044 
00046   std::string name;
00047 
00049   colvarvalue::Type type() const;
00050 
00052   colvarvalue const & value() const;
00053 
00055   colvarvalue const & velocity() const;
00056 
00065   colvarvalue const & system_force() const;
00066 
00068 
00077   cvm::real width;
00078 
00081   bool b_linear;
00082 
00085   bool b_inverse_gradients;
00086 
00089   bool b_Jacobian_force;
00090 
00093   enum task {
00096     task_gradients,
00099     task_collect_gradients,
00101     task_fdiff_velocity,
00104     task_system_force,
00108     task_extended_lagrangian,
00111     task_langevin,
00114     task_output_energy,
00118     task_Jacobian_force,
00121     task_report_Jacobian_force,
00123     task_output_value,
00125     task_output_velocity,
00127     task_output_applied_force,
00129     task_output_system_force,
00131     task_lower_boundary,
00133     task_upper_boundary,
00137     task_grid,
00140     task_lower_wall,
00143     task_upper_wall,
00145     task_runave,
00147     task_corrfunc,
00149     task_ntot
00150   };
00151 
00153   bool tasks[task_ntot];
00154 
00155 protected:
00156 
00157 
00158   /*
00159     extended:
00160     H = H_{0} + \sum_{i} 1/2*\lambda*(S_i(x(t))-s_i(t))^2 \\
00161     + \sum_{i} 1/2*m_i*(ds_i(t)/dt)^2 \\
00162     + \sum_{t'<t} W * exp (-1/2*\sum_{i} (s_i(t')-s_i(t))^2/(\delta{}s_i)^2) \\
00163     + \sum_{w} (\sum_{i}c_{w,i}s_i(t) - D_w)^M/(\Sigma_w)^M
00164 
00165     normal:
00166     H = H_{0} + \sum_{t'<t} W * exp (-1/2*\sum_{i} (S_i(x(t'))-S_i(x(t)))^2/(\delta{}S_i)^2) \\
00167     + \sum_{w} (\sum_{i}c_{w,i}S_i(t) - D_w)^M/(\Sigma_w)^M
00168 
00169     output:
00170     H = H_{0}   (only output S(x), no forces)
00171 
00172     Here:
00173     S(x(t)) = x
00174     s(t)    = xr
00175     DS = Ds = delta
00176   */
00177 
00178 
00180   colvarvalue x;
00181 
00183   colvarvalue x_reported;
00184 
00186   colvarvalue v_fdiff;
00187 
00188   inline colvarvalue fdiff_velocity (colvarvalue const &xold,
00189                                      colvarvalue const &xnew)
00190   {
00191     // using the gradient of the square distance to calculate the
00192     // velocity (non-scalar variables automatically taken into
00193     // account)
00194     cvm::real dt = cvm::dt();
00195     return ( ( (dt > 0.0) ? (1.0/dt) : 1.0 ) *
00196              0.5 * dist2_lgrad (xnew, xold) );
00197   }
00198 
00200   colvarvalue v_reported;
00201 
00202   // Options for task_extended_lagrangian
00204   colvarvalue xr;
00206   colvarvalue vr;
00208   cvm::real ext_mass;
00210   cvm::real ext_force_k;
00212   cvm::real ext_gamma;
00214   cvm::real ext_sigma;
00215   
00217   colvarvalue fr;
00218 
00220   colvarvalue fj;
00221 
00223   colvarvalue ft_reported;
00224 
00225 public:
00226 
00227 
00230   colvarvalue fb;
00231 
00235   colvarvalue f;
00236 
00239   colvarvalue ft;
00240 
00241 
00243   cvm::real period;
00244 
00247   bool b_periodic;
00248 
00249 
00252   bool expand_boundaries;
00253 
00255   colvarvalue lower_boundary;
00257   colvarvalue lower_wall;
00259   cvm::real   lower_wall_k;
00260 
00262   colvarvalue upper_boundary;
00264   colvarvalue upper_wall;
00266   cvm::real   upper_wall_k;
00267 
00269   bool periodic_boundaries() const;
00270 
00272   bool periodic_boundaries (colvarvalue const &lb, colvarvalue const &ub) const;
00273 
00274 
00276   colvar (std::string const &conf);
00277 
00279   void enable (colvar::task const &t);
00280 
00282   void disable (colvar::task const &t);
00283 
00285   ~colvar();
00286 
00287 
00290   void calc();
00291 
00293   void calc_value (colvarvalue &xn);
00294 
00296   void reset_bias_force();
00297 
00299   void add_bias_force (colvarvalue const &force);
00300 
00305   cvm::real update();
00306 
00309   void communicate_forces();
00310 
00311 
00316   cvm::real dist2 (colvarvalue const &x1,
00317                    colvarvalue const &x2) const;
00318 
00323   colvarvalue dist2_lgrad (colvarvalue const &x1,
00324                            colvarvalue const &x2) const;
00325 
00330   colvarvalue dist2_rgrad (colvarvalue const &x1,
00331                            colvarvalue const &x2) const;
00332 
00337   cvm::real compare (colvarvalue const &x1,
00338                      colvarvalue const &x2) const;
00339 
00344   void wrap (colvarvalue &x) const;
00345 
00346 
00348   void parse_analysis (std::string const &conf);
00350   void analyse();
00351 
00352 
00354   std::istream & read_traj (std::istream &is);
00356   std::ostream & write_traj (std::ostream &os);
00358   std::ostream & write_traj_label (std::ostream &os);
00359 
00360 
00362   std::istream & read_restart (std::istream &is);
00364   std::ostream & write_restart (std::ostream &os);
00365 
00366 
00367 protected:
00368 
00370   colvarvalue            x_old;
00371 
00374   std::list< std::list<colvarvalue> > acf_x_history, acf_v_history;
00377   std::list< std::list<colvarvalue> >::iterator acf_x_history_p, acf_v_history_p;
00378 
00380   std::list< std::list<colvarvalue> > x_history;
00383   std::list< std::list<colvarvalue> >::iterator x_history_p;
00384 
00387   std::string            acf_colvar_name;
00389   size_t                 acf_length;
00391   size_t                 acf_offset;
00393   size_t                 acf_stride;
00395   size_t                 acf_nframes;
00397   bool                   acf_normalize;
00399   std::vector<cvm::real> acf;
00401   std::string            acf_outfile;
00402 
00404   enum acf_type_e {
00406     acf_notset,
00408     acf_vel,
00410     acf_coor,
00413     acf_p2coor
00414   };
00415 
00417   acf_type_e             acf_type;
00418 
00420   void calc_vel_acf (std::list<colvarvalue> &v_history,
00421                      colvarvalue const      &v);
00422 
00425   void calc_coor_acf (std::list<colvarvalue> &x_history,
00426                       colvarvalue const      &x);
00427 
00430   void calc_p2coor_acf (std::list<colvarvalue> &x_history,
00431                         colvarvalue const      &x);
00432 
00434   void calc_acf();
00436   void write_acf (std::ostream &os);
00437 
00439   size_t         runave_length;
00441   size_t         runave_stride;
00443   std::ofstream  runave_os;
00445   colvarvalue    runave;
00447   cvm::real      runave_variance;
00448 
00450   void calc_runave();
00451 
00453   cvm::real kinetic_energy;
00454   cvm::real potential_energy;
00455 public:
00456 
00457 
00458   // collective variable component base class
00459   class cvc;
00460 
00461   // currently available collective variable components
00462 
00463   // scalar colvar components
00464   class distance;
00465   class distance_z;
00466   class distance_xy;
00467   class min_distance;
00468   class angle;
00469   class dihedral;
00470   class coordnum;
00471   class selfcoordnum;
00472   class h_bond;
00473   class rmsd;
00474   class logmsd;
00475   class orientation_angle;
00476   class tilt;
00477   class spin_angle;
00478   class gyration;
00479   class eigenvector;
00480   class alpha_dihedrals;
00481   class alpha_angles;
00482   class dihedPC;
00483 
00484   // non-scalar components
00485   class distance_vec;
00486   class distance_dir;
00487   class orientation;
00488 
00489 protected:
00490 
00492   std::vector<cvc *> cvcs;
00493 
00496   void build_atom_list (void);
00497 
00498 public:
00500   std::vector<int> atom_ids;
00501 
00505   std::vector<cvm::rvector> atomic_gradients;
00506 
00507   inline size_t n_components () const {
00508     return cvcs.size();
00509   }
00510 };
00511 
00512 
00513 inline colvar * cvm::colvar_p (std::string const &name)
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 }
00524 
00525 
00526 inline colvarvalue::Type colvar::type() const
00527 {
00528   return x.type();
00529 }
00530 
00531 
00532 inline colvarvalue const & colvar::value() const
00533 {
00534   return x_reported;
00535 }
00536 
00537 
00538 inline colvarvalue const & colvar::velocity() const
00539 {
00540   return v_reported;
00541 }
00542 
00543 
00544 inline colvarvalue const & colvar::system_force() const
00545 {
00546   return ft_reported;
00547 }
00548 
00549 
00550 inline void colvar::add_bias_force (colvarvalue const &force)
00551 {
00552   fb += force;
00553 }
00554 
00555 
00556 inline void colvar::reset_bias_force() {
00557   fb.reset();
00558 }
00559 
00560 #endif
00561 

Generated on Fri May 25 04:07:13 2012 for NAMD by  doxygen 1.3.9.1