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

colvarmodule.C

Go to the documentation of this file.
00001 #include "colvarmodule.h"
00002 #include "colvarparse.h"
00003 #include "colvarproxy.h"
00004 #include "colvar.h"
00005 #include "colvarbias.h"
00006 #include "colvarbias_meta.h"
00007 #include "colvarbias_abf.h"
00008 
00009 
00010 colvarmodule::colvarmodule (char const  *config_filename,
00011                             colvarproxy *proxy_in)
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 }
00136 
00137 
00138 std::istream & colvarmodule::read_restart (std::istream &is) 
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 }
00176 
00177 
00178 
00179 void colvarmodule::init_colvars (std::string const &conf)
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 }
00210 
00211 
00212 void colvarmodule::init_biases (std::string const &conf)
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 }
00298 
00299 
00300 void colvarmodule::change_configuration (std::string const &bias_name,
00301                                          std::string const &conf)
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 }
00317 
00318 cvm::real colvarmodule::energy_difference (std::string const &bias_name,
00319                                            std::string const &conf)
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 }
00337 
00338 
00339 void colvarmodule::calc() {
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 }
00514 
00515 
00516 void colvarmodule::analyze()
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 }
00544 
00545 
00546 colvarmodule::~colvarmodule()
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 }  
00568 
00569 
00570 void colvarmodule::write_output_files()
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 }
00596 
00597 
00598 
00599 bool colvarmodule::read_traj (char const *traj_filename,
00600                               size_t      traj_read_begin,
00601                               size_t      traj_read_end)
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 }
00665 
00666 
00667 std::ostream & colvarmodule::write_restart (std::ostream &os)
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 }
00691 
00692 
00693 
00694 void cvm::log (std::string const &message)
00695 {
00696   if (depth > 0)
00697     proxy->log ((std::string (2*depth, ' '))+message);
00698   else
00699     proxy->log (message);
00700 }
00701 
00702 void cvm::increase_depth()
00703 {
00704   depth++;
00705 }
00706 
00707 void cvm::decrease_depth()
00708 {
00709   if (depth) depth--;
00710 }
00711 
00712 void cvm::fatal_error (std::string const &message)
00713 {
00714   proxy->fatal_error (message);
00715 }
00716 
00717 void cvm::exit (std::string const &message)
00718 {
00719   proxy->exit (message);
00720 }
00721 
00722 
00723 
00724 // static pointers
00725 std::vector<colvar *>     colvarmodule::colvars;
00726 std::vector<colvarbias *> colvarmodule::biases;
00727 size_t                    colvarmodule::n_abf_biases = 0;
00728 size_t                    colvarmodule::n_harm_biases = 0;
00729 size_t                    colvarmodule::n_histo_biases = 0;
00730 size_t                    colvarmodule::n_meta_biases = 0;
00731 colvarproxy              *colvarmodule::proxy = NULL;
00732 
00733 
00734 // static runtime data
00735 cvm::real colvarmodule::debug_gradients_step_size = 1.0e-03;
00736 size_t    colvarmodule::it = 0;
00737 size_t    colvarmodule::it_restart = 0;
00738 size_t    colvarmodule::restart_out_freq = 0;
00739 size_t    colvarmodule::cv_traj_freq = 0;
00740 size_t    colvarmodule::depth = 0;
00741 bool      colvarmodule::b_analysis = false;
00742 cvm::real colvarmodule::rotation::crossing_threshold = 1.0E-04;
00743 
00744 
00745 // file name prefixes
00746 std::string colvarmodule::output_prefix = "";
00747 std::string colvarmodule::input_prefix = "";
00748 std::string colvarmodule::restart_in_name = "";
00749 
00750 
00751 // i/o constants
00752 size_t const colvarmodule::it_width = 12;
00753 size_t const colvarmodule::cv_prec  = 14;
00754 size_t const colvarmodule::cv_width = 21;
00755 size_t const colvarmodule::en_prec  = 14;
00756 size_t const colvarmodule::en_width = 21;
00757 std::string const colvarmodule::line_marker =
00758   "----------------------------------------------------------------------\n";
00759 
00760 
00761 
00762 
00763 
00764 
00765 std::ostream & operator << (std::ostream &os, colvarmodule::rvector const &v)
00766 {
00767   std::streamsize const w = os.width();
00768   std::streamsize const p = os.precision();
00769 
00770   os.width (2);
00771   os << "( ";
00772   os.width (w); os.precision (p);
00773   os << v.x << " , ";
00774   os.width (w); os.precision (p);
00775   os << v.y << " , ";
00776   os.width (w); os.precision (p);
00777   os << v.z << " )";
00778   return os;
00779 }
00780 
00781 
00782 std::istream & operator >> (std::istream &is, colvarmodule::rvector &v)
00783 {
00784   size_t const start_pos = is.tellg();
00785   char sep;
00786   if ( !(is >> sep) || !(sep == '(') ||
00787        !(is >> v.x) || !(is >> sep)  || !(sep == ',') ||
00788        !(is >> v.y) || !(is >> sep)  || !(sep == ',') ||
00789        !(is >> v.z) || !(is >> sep)  || !(sep == ')') ) {
00790     is.clear();
00791     is.seekg (start_pos, std::ios::beg);
00792     is.setstate (std::ios::failbit);
00793     return is;
00794   }
00795   return is;
00796 }
00797 
00798 
00799 
00800 std::ostream & operator << (std::ostream &os, colvarmodule::quaternion const &q)
00801 {
00802   std::streamsize const w = os.width();
00803   std::streamsize const p = os.precision();
00804 
00805   os.width (2);
00806   os << "( ";
00807   os.width (w); os.precision (p);
00808   os << q.q0 << " , ";
00809   os.width (w); os.precision (p);
00810   os << q.q1 << " , ";
00811   os.width (w); os.precision (p);
00812   os << q.q2 << " , ";
00813   os.width (w); os.precision (p);
00814   os << q.q3 << " )";
00815   return os;
00816 }
00817 
00818 
00819 std::istream & operator >> (std::istream &is, colvarmodule::quaternion &q)
00820 {
00821   size_t const start_pos = is.tellg();
00822 
00823   std::string euler ("");
00824 
00825   if ( (is >> euler) && (colvarparse::to_lower_cppstr (euler) ==
00826                          std::string ("euler")) ) {
00827 
00828     // parse the Euler angles
00829     
00830     char sep;
00831     cvm::real phi, theta, psi;
00832     if ( !(is >> sep)   || !(sep == '(') ||
00833          !(is >> phi)   || !(is >> sep)  || !(sep == ',') ||
00834          !(is >> theta) || !(is >> sep)  || !(sep == ',') ||
00835          !(is >> psi)   || !(is >> sep)  || !(sep == ')') ) {
00836       is.clear();
00837       is.seekg (start_pos, std::ios::beg);
00838       is.setstate (std::ios::failbit);
00839       return is;
00840     }
00841 
00842     q = colvarmodule::quaternion (phi, theta, psi);
00843 
00844   } else {
00845 
00846     // parse the quaternion components
00847 
00848     is.seekg (start_pos, std::ios::beg);      
00849     char sep;
00850     if ( !(is >> sep)  || !(sep == '(') ||
00851          !(is >> q.q0) || !(is >> sep)  || !(sep == ',') ||
00852          !(is >> q.q1) || !(is >> sep)  || !(sep == ',') ||
00853          !(is >> q.q2) || !(is >> sep)  || !(sep == ',') ||
00854          !(is >> q.q3) || !(is >> sep)  || !(sep == ')') ) {
00855       is.clear();
00856       is.seekg (start_pos, std::ios::beg);
00857       is.setstate (std::ios::failbit);
00858       return is;
00859     }
00860   }
00861 
00862   return is;
00863 }
00864 
00865 
00866 cvm::quaternion
00867 cvm::quaternion::position_derivative_inner (cvm::rvector const &pos,
00868                                             cvm::rvector const &vec) const
00869 {
00870   cvm::quaternion result (0.0, 0.0, 0.0, 0.0);
00871 
00872 
00873   result.q0 =   2.0 * pos.x * q0 * vec.x
00874                +2.0 * pos.y * q0 * vec.y
00875                +2.0 * pos.z * q0 * vec.z
00876 
00877                -2.0 * pos.y * q3 * vec.x
00878                +2.0 * pos.z * q2 * vec.x
00879 
00880                +2.0 * pos.x * q3 * vec.y
00881                -2.0 * pos.z * q1 * vec.y
00882 
00883                -2.0 * pos.x * q2 * vec.z
00884                +2.0 * pos.y * q1 * vec.z;
00885 
00886 
00887   result.q1 =  +2.0 * pos.x * q1 * vec.x
00888                -2.0 * pos.y * q1 * vec.y
00889                -2.0 * pos.z * q1 * vec.z
00890 
00891                +2.0 * pos.y * q2 * vec.x
00892                +2.0 * pos.z * q3 * vec.x
00893 
00894                +2.0 * pos.x * q2 * vec.y
00895                -2.0 * pos.z * q0 * vec.y
00896 
00897                +2.0 * pos.x * q3 * vec.z
00898                +2.0 * pos.y * q0 * vec.z;
00899 
00900 
00901   result.q2 =  -2.0 * pos.x * q2 * vec.x
00902                +2.0 * pos.y * q2 * vec.y
00903                -2.0 * pos.z * q2 * vec.z
00904 
00905                +2.0 * pos.y * q1 * vec.x
00906                +2.0 * pos.z * q0 * vec.x
00907 
00908                +2.0 * pos.x * q1 * vec.y
00909                +2.0 * pos.z * q3 * vec.y
00910 
00911                -2.0 * pos.x * q0 * vec.z
00912                +2.0 * pos.y * q3 * vec.z;
00913 
00914 
00915   result.q3 =  -2.0 * pos.x * q3 * vec.x
00916                -2.0 * pos.y * q3 * vec.y
00917                +2.0 * pos.z * q3 * vec.z
00918 
00919                -2.0 * pos.y * q0 * vec.x
00920                +2.0 * pos.z * q1 * vec.x
00921 
00922                +2.0 * pos.x * q0 * vec.y
00923                +2.0 * pos.z * q2 * vec.y
00924 
00925                +2.0 * pos.x * q1 * vec.z
00926                +2.0 * pos.y * q2 * vec.z;
00927 
00928   return result;
00929 }
00930 
00931 
00932 
00933 
00934 // Calculate the optimal rotation between two groups, and implement it
00935 // as a quaternion.  The method is the one documented in: Coutsias EA,
00936 // Seok C, Dill KA.  Using quaternions to calculate RMSD.  J Comput
00937 // Chem. 25(15):1849-57 (2004) DOI: 10.1002/jcc.20110 PubMed: 15376254
00938 
00939 
00940 
00941 
00942 
00943 void colvarmodule::rotation::build_matrix (std::vector<cvm::atom_pos> const &pos1,
00944                                            std::vector<cvm::atom_pos> const &pos2,
00945                                            matrix2d<cvm::real, 4, 4>        &S)
00946 {
00947   cvm::rmatrix C;
00948 
00949   // build the correlation matrix
00950   C.reset();
00951   for (size_t i = 0; i < pos1.size(); i++) {
00952     C.xx() += pos1[i].x * pos2[i].x;
00953     C.xy() += pos1[i].x * pos2[i].y;
00954     C.xz() += pos1[i].x * pos2[i].z;
00955     C.yx() += pos1[i].y * pos2[i].x;
00956     C.yy() += pos1[i].y * pos2[i].y;
00957     C.yz() += pos1[i].y * pos2[i].z;
00958     C.zx() += pos1[i].z * pos2[i].x;
00959     C.zy() += pos1[i].z * pos2[i].y;
00960     C.zz() += pos1[i].z * pos2[i].z;
00961   } 
00962 
00963   // build the "overlap" matrix, whose eigenvectors are stationary
00964   // points of the RMSD in the space of rotations
00965   S[0][0] =    C.xx() + C.yy() + C.zz();
00966   S[1][0] =    C.yz() - C.zy();
00967   S[0][1] = S[1][0];
00968   S[2][0] =  - C.xz() + C.zx() ;
00969   S[0][2] = S[2][0];
00970   S[3][0] =    C.xy() - C.yx();
00971   S[0][3] = S[3][0];
00972   S[1][1] =    C.xx() - C.yy() - C.zz();
00973   S[2][1] =    C.xy() + C.yx();
00974   S[1][2] = S[2][1];
00975   S[3][1] =    C.xz() + C.zx();
00976   S[1][3] = S[3][1];
00977   S[2][2] = - C.xx() + C.yy() - C.zz();
00978   S[3][2] =   C.yz() + C.zy();
00979   S[2][3] = S[3][2];
00980   S[3][3] = - C.xx() - C.yy() + C.zz();
00981 
00982   //   if (cvm::debug()) {
00983   //     for (size_t i = 0; i < 4; i++) {
00984   //       std::string line ("");
00985   //       for (size_t j = 0; j < 4; j++) {
00986   //         line += std::string (" S["+cvm::to_str (i)+
00987   //                              "]["+cvm::to_str (j)+"] ="+cvm::to_str (S[i][j]));
00988   //       }
00989   //       cvm::log (line+"\n");
00990   //     }
00991   //   }
00992 }
00993 
00994 
00995 void colvarmodule::rotation::diagonalize_matrix (matrix2d<cvm::real, 4, 4> &S,
00996                                                  cvm::real                  S_eigval[4],
00997                                                  matrix2d<cvm::real, 4, 4> &S_eigvec)
00998 {
00999   // diagonalize
01000   int jac_nrot = 0;
01001   jacobi (S, 4, S_eigval, S_eigvec, &jac_nrot);
01002   eigsrt (S_eigval, S_eigvec, 4);
01003   // jacobi saves eigenvectors by columns
01004   transpose (S_eigvec, 4);
01005 
01006   // normalize eigenvectors
01007   for (size_t ie = 0; ie < 4; ie++) {
01008     cvm::real norm2 = 0.0;
01009     for (size_t i = 0; i < 4; i++) norm2 += std::pow (S_eigvec[ie][i], int (2));
01010     cvm::real const norm = std::sqrt (norm2);
01011     for (size_t i = 0; i < 4; i++) S_eigvec[ie][i] /= norm;
01012   }
01013 }
01014 
01015 
01016 // Calculate the rotation, plus its derivatives
01017 
01018 void colvarmodule::rotation::calc_optimal_rotation
01019 (std::vector<cvm::atom_pos> const &pos1,
01020  std::vector<cvm::atom_pos> const &pos2)
01021 {
01022   matrix2d<cvm::real, 4, 4> S;
01023   matrix2d<cvm::real, 4, 4> S_backup;
01024   cvm::real                 S_eigval[4];
01025   matrix2d<cvm::real, 4, 4> S_eigvec;
01026 
01027 //   if (cvm::debug()) {
01028 //     cvm::atom_pos cog1 (0.0, 0.0, 0.0);
01029 //     for (size_t i = 0; i < pos1.size(); i++) {
01030 //       cog1 += pos1[i];
01031 //     }
01032 //     cog1 /= cvm::real (pos1.size());
01033 //     cvm::atom_pos cog2 (0.0, 0.0, 0.0);
01034 //     for (size_t i = 0; i < pos2.size(); i++) {
01035 //       cog2 += pos2[i];
01036 //     }
01037 //     cog2 /= cvm::real (pos1.size());
01038 //     cvm::log ("calc_optimal_rotation: centers of geometry are: "+
01039 //               cvm::to_str (cog1, cvm::cv_width, cvm::cv_prec)+
01040 //               " and "+cvm::to_str (cog2, cvm::cv_width, cvm::cv_prec)+".\n");
01041 //   }
01042 
01043   build_matrix (pos1, pos2, S);
01044   S_backup = S;
01045 
01046   if (cvm::debug()) {
01047     if (b_debug_gradients) {
01048       cvm::log ("S     = "+cvm::to_str (cvm::to_str (S_backup), cvm::cv_width, cvm::cv_prec)+"\n");
01049     }
01050   }
01051 
01052   diagonalize_matrix (S, S_eigval, S_eigvec);
01053 
01054   // eigenvalues and eigenvectors
01055   cvm::real const &L0 = S_eigval[0];
01056   cvm::real const &L1 = S_eigval[1];
01057   cvm::real const &L2 = S_eigval[2];
01058   cvm::real const &L3 = S_eigval[3];
01059   cvm::real const *Q0 = S_eigvec[0];
01060   cvm::real const *Q1 = S_eigvec[1];
01061   cvm::real const *Q2 = S_eigvec[2];
01062   cvm::real const *Q3 = S_eigvec[3];
01063 
01064   lambda = L0;
01065   q = cvm::quaternion (Q0);
01066 
01067   if (q_old.norm2() > 0.0) {
01068     q.match (q_old);
01069     if (q_old.inner (q) < (1.0 - crossing_threshold)) {
01070       cvm::log ("Warning: discontinuous rotation!\n");
01071     }
01072   }
01073   q_old = q;
01074                      
01075   if (cvm::debug()) {
01076     if (b_debug_gradients) {
01077       cvm::log ("L0 = "+cvm::to_str (L0, cvm::cv_width, cvm::cv_prec)+
01078                 ", Q0 = "+cvm::to_str (cvm::quaternion (Q0), cvm::cv_width, cvm::cv_prec)+
01079                 ", Q0*Q0 = "+cvm::to_str (cvm::quaternion (Q0).inner (cvm::quaternion (Q0)), cvm::cv_width, cvm::cv_prec)+
01080                 "\n");
01081       cvm::log ("L1 = "+cvm::to_str (L1, cvm::cv_width, cvm::cv_prec)+
01082                 ", Q1 = "+cvm::to_str (cvm::quaternion (Q1), cvm::cv_width, cvm::cv_prec)+
01083                 ", Q0*Q1 = "+cvm::to_str (cvm::quaternion (Q0).inner (cvm::quaternion (Q1)), cvm::cv_width, cvm::cv_prec)+
01084                 "\n");
01085       cvm::log ("L2 = "+cvm::to_str (L2, cvm::cv_width, cvm::cv_prec)+
01086                 ", Q2 = "+cvm::to_str (cvm::quaternion (Q2), cvm::cv_width, cvm::cv_prec)+
01087                 ", Q0*Q2 = "+cvm::to_str (cvm::quaternion (Q0).inner (cvm::quaternion (Q2)), cvm::cv_width, cvm::cv_prec)+
01088                 "\n");
01089       cvm::log ("L3 = "+cvm::to_str (L3, cvm::cv_width, cvm::cv_prec)+
01090                 ", Q3 = "+cvm::to_str (cvm::quaternion (Q3), cvm::cv_width, cvm::cv_prec)+
01091                 ", Q0*Q3 = "+cvm::to_str (cvm::quaternion (Q0).inner (cvm::quaternion (Q3)), cvm::cv_width, cvm::cv_prec)+
01092                 "\n");
01093     }
01094   }
01095 
01096   // calculate derivatives of L0 and Q0 with respect to each atom in
01097   // either group; note: if dS_1 is a null vector, nothing will be
01098   // calculated
01099   for (size_t ia = 0; ia < dS_1.size(); ia++) {
01100 
01101     cvm::real const &a2x = pos2[ia].x;
01102     cvm::real const &a2y = pos2[ia].y;
01103     cvm::real const &a2z = pos2[ia].z;
01104 
01105     matrix2d<cvm::rvector, 4, 4>    &ds_1 =  dS_1[ia];
01106 
01107     // derivative of the S matrix
01108     ds_1.reset();
01109     ds_1[0][0] = cvm::rvector ( a2x,  a2y,  a2z);
01110     ds_1[1][0] = cvm::rvector ( 0.0,  a2z, -a2y);
01111     ds_1[0][1] = ds_1[1][0];
01112     ds_1[2][0] = cvm::rvector (-a2z,  0.0,  a2x);
01113     ds_1[0][2] = ds_1[2][0];
01114     ds_1[3][0] = cvm::rvector ( a2y, -a2x,  0.0);
01115     ds_1[0][3] = ds_1[3][0];
01116     ds_1[1][1] = cvm::rvector ( a2x, -a2y, -a2z);
01117     ds_1[2][1] = cvm::rvector ( a2y,  a2x,  0.0);
01118     ds_1[1][2] = ds_1[2][1];
01119     ds_1[3][1] = cvm::rvector ( a2z,  0.0,  a2x);
01120     ds_1[1][3] = ds_1[3][1];
01121     ds_1[2][2] = cvm::rvector (-a2x,  a2y, -a2z);
01122     ds_1[3][2] = cvm::rvector ( 0.0,  a2z,  a2y);
01123     ds_1[2][3] = ds_1[3][2];
01124     ds_1[3][3] = cvm::rvector (-a2x, -a2y,  a2z);
01125 
01126     cvm::rvector              &dl0_1 = dL0_1[ia];
01127     vector1d<cvm::rvector, 4> &dq0_1 = dQ0_1[ia];
01128 
01129     // matrix multiplications; derivatives of L_0 and Q_0 are
01130     // calculated using Hellmann-Feynman theorem (i.e. exploiting the
01131     // fact that the eigenvectors Q_i form an orthonormal basis)
01132 
01133     dl0_1.reset();
01134     for (size_t i = 0; i < 4; i++) {
01135       for (size_t j = 0; j < 4; j++) {
01136         dl0_1 += Q0[i] * ds_1[i][j] * Q0[j];
01137       }
01138     }
01139 
01140     dq0_1.reset();
01141     for (size_t p = 0; p < 4; p++) {
01142       for (size_t i = 0; i < 4; i++) {
01143         for (size_t j = 0; j < 4; j++) {
01144           dq0_1[p] +=
01145             (Q1[i] * ds_1[i][j] * Q0[j]) / (L0-L1) * Q1[p] + 
01146             (Q2[i] * ds_1[i][j] * Q0[j]) / (L0-L2) * Q2[p] + 
01147             (Q3[i] * ds_1[i][j] * Q0[j]) / (L0-L3) * Q3[p];
01148         }
01149       }
01150     }
01151   }
01152 
01153   // do the same for the second group
01154   for (size_t ia = 0; ia < dS_2.size(); ia++) {
01155 
01156     cvm::real const &a1x = pos1[ia].x;
01157     cvm::real const &a1y = pos1[ia].y;
01158     cvm::real const &a1z = pos1[ia].z;
01159 
01160     matrix2d<cvm::rvector, 4, 4> &ds_2 =  dS_2[ia];
01161 
01162     ds_2.reset();
01163     ds_2[0][0] = cvm::rvector ( a1x,  a1y,  a1z);
01164     ds_2[1][0] = cvm::rvector ( 0.0, -a1z,  a1y);
01165     ds_2[0][1] = ds_2[1][0];
01166     ds_2[2][0] = cvm::rvector ( a1z,  0.0, -a1x);
01167     ds_2[0][2] = ds_2[2][0];
01168     ds_2[3][0] = cvm::rvector (-a1y,  a1x,  0.0);
01169     ds_2[0][3] = ds_2[3][0];
01170     ds_2[1][1] = cvm::rvector ( a1x, -a1y, -a1z);
01171     ds_2[2][1] = cvm::rvector ( a1y,  a1x,  0.0);
01172     ds_2[1][2] = ds_2[2][1];
01173     ds_2[3][1] = cvm::rvector ( a1z,  0.0,  a1x);
01174     ds_2[1][3] = ds_2[3][1];
01175     ds_2[2][2] = cvm::rvector (-a1x,  a1y, -a1z);
01176     ds_2[3][2] = cvm::rvector ( 0.0,  a1z,  a1y);
01177     ds_2[2][3] = ds_2[3][2];
01178     ds_2[3][3] = cvm::rvector (-a1x, -a1y,  a1z);
01179 
01180     cvm::rvector              &dl0_2 = dL0_2[ia];
01181     vector1d<cvm::rvector, 4> &dq0_2 = dQ0_2[ia];
01182 
01183     dl0_2.reset();
01184     for (size_t i = 0; i < 4; i++) {
01185       for (size_t j = 0; j < 4; j++) {
01186         dl0_2 += Q0[i] * ds_2[i][j] * Q0[j];
01187       }
01188     }
01189 
01190     dq0_2.reset();
01191     for (size_t p = 0; p < 4; p++) {
01192       for (size_t i = 0; i < 4; i++) {
01193         for (size_t j = 0; j < 4; j++) {
01194           dq0_2[p] +=
01195             (Q1[i] * ds_2[i][j] * Q0[j]) / (L0-L1) * Q1[p] + 
01196             (Q2[i] * ds_2[i][j] * Q0[j]) / (L0-L2) * Q2[p] + 
01197             (Q3[i] * ds_2[i][j] * Q0[j]) / (L0-L3) * Q3[p];
01198         }
01199       }
01200     }
01201 
01202     if (cvm::debug()) {
01203       
01204       if (b_debug_gradients) {
01205 
01206       matrix2d<cvm::real, 4, 4> S_new;
01207       cvm::real                 S_new_eigval[4];
01208       matrix2d<cvm::real, 4, 4> S_new_eigvec;
01209 
01210         // make an infitesimal move along each cartesian coordinate of
01211         // this atom, and solve again the eigenvector problem
01212         for (size_t comp = 0; comp < 3; comp++) {
01213 
01214           S_new = S_backup;
01215           // diagonalize the new overlap matrix
01216           for (size_t i = 0; i < 4; i++) {
01217             for (size_t j = 0; j < 4; j++) {
01218               S_new[i][j] += 
01219                 colvarmodule::debug_gradients_step_size * ds_2[i][j][comp];
01220             }
01221           }
01222 
01223 //           cvm::log ("S_new = "+cvm::to_str (cvm::to_str (S_new), cvm::cv_width, cvm::cv_prec)+"\n");
01224 
01225           diagonalize_matrix (S_new, S_new_eigval, S_new_eigvec);
01226 
01227           cvm::real const &L0_new = S_new_eigval[0];
01228           cvm::real const *Q0_new = S_new_eigvec[0];
01229 
01230           cvm::real const DL0 = (dl0_2[comp]) * colvarmodule::debug_gradients_step_size;
01231           cvm::quaternion const q0 (Q0);
01232           cvm::quaternion const DQ0 (dq0_2[0][comp] * colvarmodule::debug_gradients_step_size,
01233                                      dq0_2[1][comp] * colvarmodule::debug_gradients_step_size,
01234                                      dq0_2[2][comp] * colvarmodule::debug_gradients_step_size,
01235                                      dq0_2[3][comp] * colvarmodule::debug_gradients_step_size);
01236 
01237           cvm::log (  "|(l_0+dl_0) - l_0^new|/l_0 = "+
01238                       cvm::to_str (std::fabs (L0+DL0 - L0_new)/L0, cvm::cv_width, cvm::cv_prec)+
01239                       ", |(q_0+dq_0) - q_0^new| = "+
01240                       cvm::to_str ((Q0+DQ0 - Q0_new).norm(), cvm::cv_width, cvm::cv_prec)+
01241                       "\n");
01242         }
01243       }
01244     }
01245   }
01246 }
01247 
01248 
01249 
01250 // Numerical Recipes routine for diagonalization
01251 
01252 #define ROTATE(a,i,j,k,l) g=a[i][j];h=a[k][l];a[i][j]=g-s*(h+g*tau);    \
01253   a[k][l]=h+s*(g-h*tau);
01254 void jacobi(cvm::real **a, int n, cvm::real d[], cvm::real **v, int *nrot)
01255 {
01256   int j,iq,ip,i;
01257   cvm::real tresh,theta,tau,t,sm,s,h,g,c;
01258 
01259   std::vector<cvm::real> b (n, 0.0);
01260   std::vector<cvm::real> z (n, 0.0);
01261 
01262   for (ip=0;ip<n;ip++) {
01263     for (iq=0;iq<n;iq++) v[ip][iq]=0.0;
01264     v[ip][ip]=1.0;
01265   }
01266   for (ip=0;ip<n;ip++) {
01267     b[ip]=d[ip]=a[ip][ip];
01268     z[ip]=0.0;
01269   }
01270   *nrot=0;
01271   for (i=0;i<=50;i++) {
01272     sm=0.0;
01273     for (ip=0;ip<n-1;ip++) {
01274       for (iq=ip+1;iq<n;iq++)
01275         sm += std::fabs(a[ip][iq]);
01276     }
01277     if (sm == 0.0) {
01278       return;
01279     }
01280     if (i < 4)
01281       tresh=0.2*sm/(n*n);
01282     else
01283       tresh=0.0;
01284     for (ip=0;ip<n-1;ip++) {
01285       for (iq=ip+1;iq<n;iq++) {
01286         g=100.0*std::fabs(a[ip][iq]);
01287         if (i > 4 && (cvm::real)(std::fabs(d[ip])+g) == (cvm::real)std::fabs(d[ip])
01288             && (cvm::real)(std::fabs(d[iq])+g) == (cvm::real)std::fabs(d[iq]))
01289           a[ip][iq]=0.0;
01290         else if (std::fabs(a[ip][iq]) > tresh) {
01291           h=d[iq]-d[ip];
01292           if ((cvm::real)(std::fabs(h)+g) == (cvm::real)std::fabs(h))
01293             t=(a[ip][iq])/h;
01294           else {
01295             theta=0.5*h/(a[ip][iq]);
01296             t=1.0/(std::fabs(theta)+std::sqrt(1.0+theta*theta));
01297             if (theta < 0.0) t = -t;
01298           }
01299           c=1.0/std::sqrt(1+t*t);
01300           s=t*c;
01301           tau=s/(1.0+c);
01302           h=t*a[ip][iq];
01303           z[ip] -= h;
01304           z[iq] += h;
01305           d[ip] -= h;
01306           d[iq] += h;
01307           a[ip][iq]=0.0;
01308           for (j=0;j<=ip-1;j++) {
01309             ROTATE(a,j,ip,j,iq)
01310               }
01311           for (j=ip+1;j<=iq-1;j++) {
01312             ROTATE(a,ip,j,j,iq)
01313               }
01314           for (j=iq+1;j<n;j++) {
01315             ROTATE(a,ip,j,iq,j)
01316               }
01317           for (j=0;j<n;j++) {
01318             ROTATE(v,j,ip,j,iq)
01319               }
01320           ++(*nrot);
01321         }
01322       }
01323     }
01324     for (ip=0;ip<n;ip++) {
01325       b[ip] += z[ip];
01326       d[ip]=b[ip];
01327       z[ip]=0.0;
01328     }
01329   }
01330   cvm::fatal_error ("Too many iterations in routine jacobi.\n");
01331 }
01332 #undef ROTATE
01333 
01334 
01335 void eigsrt(cvm::real d[], cvm::real **v, int n)
01336 {
01337   int k,j,i;
01338   cvm::real p;
01339 
01340   for (i=0;i<n;i++) {
01341     p=d[k=i];
01342     for (j=i+1;j<n;j++)
01343       if (d[j] >= p) p=d[k=j];
01344     if (k != i) {
01345       d[k]=d[i];
01346       d[i]=p;
01347       for (j=0;j<n;j++) {
01348         p=v[j][i];
01349         v[j][i]=v[j][k];
01350         v[j][k]=p;
01351       }
01352     }
01353   }
01354 }
01355 
01356 
01357 void transpose(cvm::real **v, int n)
01358 {
01359   cvm::real p;
01360   for (int i=0;i<n;i++) {
01361     for (int j=i+1;j<n;j++) {
01362       p=v[i][j];
01363       v[i][j]=v[j][i];
01364       v[j][i]=p;
01365     }
01366   }
01367 }
01368 

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