00001 #include <iostream>
00002 #include <sstream>
00003 #include <fstream>
00004 #include <iomanip>
00005 #include <cmath>
00006 #include <algorithm>
00007
00008
00009 #if defined (WIN32) && !defined(__CYGWIN__)
00010 #include <direct.h>
00011 #define CHDIR ::_chdir
00012 #define GETCWD ::_getcwd
00013 #define PATHSEP "\\"
00014 #else
00015 #include <unistd.h>
00016 #define CHDIR ::chdir
00017 #define GETCWD ::getcwd
00018 #define PATHSEP "/"
00019 #endif
00020
00021
00022 #include "colvar.h"
00023 #include "colvarbias_meta.h"
00024
00025
00026 colvarbias_meta::colvarbias_meta()
00027 : colvarbias(),
00028 state_file_step (0),
00029 new_hills_begin (hills.end())
00030 {
00031 }
00032
00033
00034 colvarbias_meta::colvarbias_meta (std::string const &conf, char const *key)
00035 : colvarbias (conf, key),
00036 state_file_step (0),
00037 new_hills_begin (hills.end())
00038 {
00039 if (cvm::n_abf_biases > 0)
00040 cvm::log ("Warning: ABF and metadynamics running at the "
00041 "same time can lead to inconsistent results.\n");
00042
00043 get_keyval (conf, "hillWeight", hill_weight, 0.01);
00044 if (hill_weight == 0.0)
00045 cvm::log ("Warning: hillWeight has been set to zero, "
00046 "this bias will have no effect.\n");
00047
00048 get_keyval (conf, "newHillFrequency", new_hill_freq, 1000);
00049
00050 get_keyval (conf, "hillWidth", hill_width, std::sqrt (2.0 * PI) / 2.0);
00051
00052 get_keyval (conf, "useGrids", use_grids, true);
00053
00054 if (use_grids) {
00055 get_keyval (conf, "gridsUpdateFrequency", grids_freq, new_hill_freq);
00056 get_keyval (conf, "rebinGrids", rebin_grids, false);
00057
00058 expand_grids = false;
00059 for (size_t i = 0; i < colvars.size(); i++) {
00060 if (colvars[i]->expand_boundaries) {
00061 expand_grids = true;
00062 cvm::log ("Metadynamics bias \""+this->name+"\""+
00063 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00064 ": Will expand grids when the colvar \""+
00065 colvars[i]->name+"\" approaches its boundaries.\n");
00066 }
00067 }
00068
00069 get_keyval (conf, "keepHills", keep_hills, false);
00070 get_keyval (conf, "dumpFreeEnergyFile", dump_fes, true);
00071 get_keyval (conf, "saveFreeEnergyFile", dump_fes_save, false);
00072
00073 for (size_t i = 0; i < colvars.size(); i++) {
00074 colvars[i]->enable (colvar::task_grid);
00075 }
00076
00077 hills_energy = new colvar_grid_scalar (colvars);
00078 hills_energy_gradients = new colvar_grid_gradient (colvars);
00079 } else {
00080 rebin_grids = false;
00081 keep_hills = false;
00082 dump_fes = false;
00083 dump_fes_save = false;
00084 dump_replica_fes = false;
00085 }
00086
00087 {
00088 bool b_replicas = false;
00089 get_keyval (conf, "multipleReplicas", b_replicas, false);
00090 if (b_replicas)
00091 comm = multiple_replicas;
00092 else
00093 comm = single_replica;
00094 }
00095
00096 if (comm != single_replica) {
00097
00098 if (expand_grids)
00099 cvm::fatal_error ("Error: expandBoundaries is not supported when "
00100 "using more than one replicas; please allocate "
00101 "wide enough boundaries for each colvar"
00102 "ahead of time.\n");
00103
00104 if (get_keyval (conf, "dumpPartialFreeEnergyFile", dump_replica_fes, false)) {
00105 if (dump_replica_fes && (! dump_fes)) {
00106 cvm::log ("Enabling \"dumpFreeEnergyFile\".\n");
00107 }
00108 }
00109
00110 get_keyval (conf, "replicaID", replica_id, std::string (""));
00111 if (!replica_id.size())
00112 cvm::fatal_error ("Error: replicaID must be defined "
00113 "when using more than one replica.\n");
00114
00115 get_keyval (conf, "replicasRegistry",
00116 replicas_registry_file,
00117 (this->name+".replicas.registry.txt"));
00118
00119 get_keyval (conf, "replicaUpdateFrequency",
00120 replica_update_freq, new_hill_freq);
00121
00122 if (keep_hills)
00123 cvm::log ("Warning: in metadynamics bias \""+this->name+"\""+
00124 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00125 ": keepHills with more than one replica can lead to a very "
00126 "large amount input/output and slow down your calculations. "
00127 "Please consider disabling it.\n");
00128
00129
00130 {
00131
00132 char *pwd = new char[3001];
00133 if (GETCWD (pwd, 3000) == NULL)
00134 cvm::fatal_error ("Error: cannot get the path of the current working directory.\n");
00135 replica_list_file =
00136 (std::string (pwd)+std::string (PATHSEP)+
00137 this->name+"."+replica_id+".files.txt");
00138
00139
00140
00141 replica_hills_file =
00142 (std::string (pwd)+std::string (PATHSEP)+
00143 cvm::output_prefix+".colvars."+this->name+"."+replica_id+".hills");
00144 replica_state_file =
00145 (std::string (pwd)+std::string (PATHSEP)+
00146 cvm::output_prefix+".colvars."+this->name+"."+replica_id+".state");
00147 delete pwd;
00148 }
00149
00150
00151
00152
00153 bool registered_replica = false;
00154 std::ifstream reg_is (replicas_registry_file.c_str());
00155 if (reg_is.good()) {
00156 std::string existing_replica ("");
00157 std::string existing_replica_file ("");
00158 while ((reg_is >> existing_replica) && existing_replica.size() &&
00159 (reg_is >> existing_replica_file) && existing_replica_file.size()) {
00160 if (existing_replica == replica_id) {
00161
00162 replica_list_file = existing_replica_file;
00163 reg_is.close();
00164 registered_replica = true;
00165 break;
00166 }
00167 }
00168 reg_is.close();
00169 }
00170
00171
00172
00173
00174
00175
00176 replica_hills_os.open (replica_hills_file.c_str());
00177 if (!replica_hills_os.good())
00178 cvm::fatal_error ("Error: in opening file \""+
00179 replica_hills_file+"\" for writing.\n");
00180 replica_hills_os.setf (std::ios::scientific, std::ios::floatfield);
00181
00182
00183 write_replica_state_file();
00184
00185 for (size_t ir = 0; ir < replicas.size(); ir++) {
00186 (replicas[ir])->replica_state_file_in_sync = false;
00187 }
00188
00189
00190
00191 std::ofstream list_os (replica_list_file.c_str(),
00192 (use_grids ? std::ios::trunc : std::ios::app));
00193 if (! list_os.good())
00194 cvm::fatal_error ("Error: in opening file \""+
00195 replica_list_file+"\" for writing.\n");
00196 list_os << "stateFile " << replica_state_file << "\n";
00197 list_os << "hillsFile " << replica_hills_file << "\n";
00198 list_os.close();
00199
00200
00201 if (! registered_replica) {
00202 std::ofstream reg_os (replicas_registry_file.c_str(), std::ios::app);
00203 if (! reg_os.good())
00204 cvm::fatal_error ("Error: in opening file \""+
00205 replicas_registry_file+"\" for writing.\n");
00206 reg_os << replica_id << " " << replica_list_file << "\n";
00207 reg_os.close();
00208 }
00209 }
00210
00211 get_keyval (conf, "writeHillsTrajectory", b_hills_traj, false);
00212 if (b_hills_traj) {
00213 std::string const traj_file_name (cvm::output_prefix+
00214 ".colvars."+this->name+
00215 ( (comm != single_replica) ?
00216 ("."+replica_id) :
00217 ("") )+
00218 ".hills.traj");
00219 hills_traj_os.open (traj_file_name.c_str());
00220 if (!hills_traj_os.good())
00221 cvm::fatal_error ("Error: in opening hills output file \"" +
00222 traj_file_name + "\".\n");
00223 }
00224
00225
00226 get_keyval (conf, "wellTempered", well_tempered, false);
00227 get_keyval (conf, "biasTemperature", bias_temperature, -1.0);
00228 if ((bias_temperature == -1.0) && well_tempered) {
00229 cvm::fatal_error ("Error: biasTemperature is not set.\n");
00230 }
00231 if (well_tempered) {
00232 cvm::log("Well-tempered metadynamics is used.\n");
00233 cvm::log("The bias temperature is "+cvm::to_str(bias_temperature)+".\n");
00234 }
00235
00236 if (cvm::debug())
00237 cvm::log ("Done initializing the metadynamics bias \""+this->name+"\""+
00238 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+".\n");
00239
00240 save_delimiters = false;
00241 }
00242
00243
00244 colvarbias_meta::~colvarbias_meta()
00245 {
00246 if (hills_energy) {
00247 delete hills_energy;
00248 hills_energy = NULL;
00249 }
00250
00251 if (hills_energy_gradients) {
00252 delete hills_energy_gradients;
00253 hills_energy_gradients = NULL;
00254 }
00255
00256 if (replica_hills_os.good())
00257 replica_hills_os.close();
00258
00259 if (hills_traj_os.good())
00260 hills_traj_os.close();
00261 }
00262
00263
00264
00265
00266
00267
00268
00269 std::list<colvarbias_meta::hill>::const_iterator
00270 colvarbias_meta::create_hill (colvarbias_meta::hill const &h)
00271 {
00272 hill_iter const hills_end = hills.end();
00273 hills.push_back (h);
00274 if (new_hills_begin == hills_end) {
00275
00276 new_hills_begin = hills.end();
00277 new_hills_begin--;
00278 }
00279
00280 if (use_grids) {
00281
00282
00283
00284
00285 cvm::real const min_dist = hills_energy->bin_distance_from_boundaries (h.centers);
00286 if (min_dist < (3.0 * std::floor (hill_width)) + 1.0) {
00287 hills_off_grid.push_back (h);
00288 }
00289 }
00290
00291
00292 if (hills_traj_os.good()) {
00293 hills_traj_os << (hills.back()).output_traj();
00294 if (cvm::debug()) {
00295 hills_traj_os.flush();
00296 }
00297 }
00298
00299 has_data = true;
00300 return hills.end();
00301 }
00302
00303
00304 std::list<colvarbias_meta::hill>::const_iterator
00305 colvarbias_meta::delete_hill (hill_iter &h)
00306 {
00307 if (cvm::debug()) {
00308 cvm::log ("Deleting hill from the metadynamics bias \""+this->name+"\""+
00309 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00310 ", with step number "+
00311 cvm::to_str (h->it)+(h->replica.size() ?
00312 ", replica id \""+h->replica :
00313 "")+".\n");
00314 }
00315
00316 if (use_grids && hills_off_grid.size()) {
00317 for (hill_iter hoff = hills_off_grid.begin();
00318 hoff != hills_off_grid.end(); hoff++) {
00319 if (*h == *hoff) {
00320 hills_off_grid.erase (hoff);
00321 break;
00322 }
00323 }
00324 }
00325
00326 if (hills_traj_os.good()) {
00327
00328 hills_traj_os << "# DELETED this hill: "
00329 << (hills.back()).output_traj()
00330 << "\n";
00331 if (cvm::debug())
00332 hills_traj_os.flush();
00333 }
00334
00335 return hills.erase (h);
00336 }
00337
00338
00339 cvm::real colvarbias_meta::update()
00340 {
00341 if (cvm::debug())
00342 cvm::log ("Updating the metadynamics bias \""+this->name+"\""+
00343 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+".\n");
00344
00345 if (use_grids) {
00346
00347 std::vector<int> curr_bin = hills_energy->get_colvars_index();
00348
00349 if (expand_grids) {
00350
00351
00352 if (cvm::debug())
00353 cvm::log ("Metadynamics bias \""+this->name+"\""+
00354 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00355 ": current coordinates on the grid: "+
00356 cvm::to_str (curr_bin)+".\n");
00357
00358 bool changed_grids = false;
00359 int const min_buffer =
00360 (3 * (size_t) std::floor (hill_width)) + 1;
00361
00362 std::vector<int> new_sizes (hills_energy->sizes());
00363 std::vector<colvarvalue> new_lower_boundaries (hills_energy->lower_boundaries);
00364 std::vector<colvarvalue> new_upper_boundaries (hills_energy->upper_boundaries);
00365
00366 for (size_t i = 0; i < colvars.size(); i++) {
00367
00368 if (! colvars[i]->expand_boundaries)
00369 continue;
00370
00371 cvm::real &new_lb = new_lower_boundaries[i].real_value;
00372 cvm::real &new_ub = new_upper_boundaries[i].real_value;
00373 int &new_size = new_sizes[i];
00374 bool changed_lb = false, changed_ub = false;
00375
00376 if (curr_bin[i] < min_buffer) {
00377 int const extra_points = (min_buffer - curr_bin[i]);
00378 new_lb -= extra_points * colvars[i]->width;
00379 new_size += extra_points;
00380
00381
00382 curr_bin[i] += extra_points;
00383
00384 changed_lb = true;
00385 cvm::log ("Metadynamics bias \""+this->name+"\""+
00386 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00387 ": new lower boundary for colvar \""+
00388 colvars[i]->name+"\", at "+
00389 cvm::to_str (new_lower_boundaries[i])+".\n");
00390 }
00391
00392 if (curr_bin[i] > new_size - min_buffer - 1) {
00393 int const extra_points = (curr_bin[i] - (new_size - 1) + min_buffer);
00394 new_ub += extra_points * colvars[i]->width;
00395 new_size += extra_points;
00396
00397 changed_ub = true;
00398 cvm::log ("Metadynamics bias \""+this->name+"\""+
00399 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00400 ": new upper boundary for colvar \""+
00401 colvars[i]->name+"\", at "+
00402 cvm::to_str (new_upper_boundaries[i])+".\n");
00403 }
00404
00405 if (changed_lb || changed_ub)
00406 changed_grids = true;
00407 }
00408
00409 if (changed_grids) {
00410
00411
00412
00413 colvar_grid_scalar *new_hills_energy =
00414 new colvar_grid_scalar (*hills_energy);
00415 colvar_grid_gradient *new_hills_energy_gradients =
00416 new colvar_grid_gradient (*hills_energy_gradients);
00417
00418
00419
00420 new_hills_energy->lower_boundaries = new_lower_boundaries;
00421 new_hills_energy->upper_boundaries = new_upper_boundaries;
00422 new_hills_energy->create (new_sizes, 0.0, 1);
00423
00424 new_hills_energy_gradients->lower_boundaries = new_lower_boundaries;
00425 new_hills_energy_gradients->upper_boundaries = new_upper_boundaries;
00426 new_hills_energy_gradients->create (new_sizes, 0.0, colvars.size());
00427
00428 new_hills_energy->map_grid (*hills_energy);
00429 new_hills_energy_gradients->map_grid (*hills_energy_gradients);
00430
00431 delete hills_energy;
00432 delete hills_energy_gradients;
00433 hills_energy = new_hills_energy;
00434 hills_energy_gradients = new_hills_energy_gradients;
00435
00436 curr_bin = hills_energy->get_colvars_index();
00437 if (cvm::debug())
00438 cvm::log ("Coordinates on the new grid: "+
00439 cvm::to_str (curr_bin)+".\n");
00440 }
00441 }
00442 }
00443
00444
00445 if ((cvm::step_absolute() % new_hill_freq) == 0) {
00446
00447 if (cvm::debug())
00448 cvm::log ("Metadynamics bias \""+this->name+"\""+
00449 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00450 ": adding a new hill at step "+cvm::to_str (cvm::step_absolute())+".\n");
00451
00452 switch (comm) {
00453
00454 case single_replica:
00455 if (well_tempered) {
00456 std::vector<int> curr_bin = hills_energy->get_colvars_index();
00457 cvm::real hills_energy_sum_here = hills_energy->value(curr_bin);
00458 cvm::real exp_weight = std::exp(-hills_energy_sum_here/(bias_temperature*cvm::boltzmann()));
00459 create_hill (hill ((hill_weight*exp_weight), colvars, hill_width));
00460 cvm::log ("Hill_height at step "+cvm::to_str (cvm::step_absolute())+" is "+
00461 cvm::to_str(hill_weight*exp_weight)+".\n");
00462
00463 } else {
00464 create_hill (hill (hill_weight, colvars, hill_width));
00465 }
00466 break;
00467
00468 case multiple_replicas:
00469 if (well_tempered) {
00470 std::vector<int> curr_bin = hills_energy->get_colvars_index();
00471 cvm::real hills_energy_sum_here = hills_energy->value(curr_bin);
00472 cvm::real exp_weight = std::exp(-hills_energy_sum_here/(bias_temperature*cvm::boltzmann()));
00473 create_hill (hill ((hill_weight*exp_weight), colvars, hill_width, replica_id));
00474 cvm::log ("Hill_height at step "+cvm::to_str (cvm::step_absolute())+" is "+
00475 cvm::to_str(hill_weight*exp_weight)+".\n");
00476
00477 } else {
00478 create_hill (hill (hill_weight, colvars, hill_width, replica_id));
00479 }
00480
00481 if (replica_hills_os.good()) {
00482 replica_hills_os << hills.back();
00483 } else {
00484 cvm::fatal_error ("Error: in metadynamics bias \""+this->name+"\""+
00485 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00486 " while writing hills for the other replicas.\n");
00487 }
00488 break;
00489 }
00490 }
00491
00492
00493 if (comm != single_replica) {
00494
00495
00496 if ((cvm::step_absolute() % replica_update_freq) == 0) {
00497 update_replicas_registry();
00498
00499 replica_hills_os.flush();
00500
00501 read_replica_files();
00502 }
00503 }
00504
00505
00506 bias_energy = 0.0;
00507 for (size_t ir = 0; ir < colvars.size(); ir++) {
00508 colvar_forces[ir].reset();
00509 }
00510 if (comm == multiple_replicas)
00511 for (size_t ir = 0; ir < replicas.size(); ir++) {
00512 replicas[ir]->bias_energy = 0.0;
00513 for (size_t ic = 0; ic < colvars.size(); ic++) {
00514 replicas[ir]->colvar_forces[ic].reset();
00515 }
00516 }
00517
00518 if (use_grids) {
00519
00520
00521
00522 if ((cvm::step_absolute() % grids_freq) == 0) {
00523
00524 project_hills (new_hills_begin, hills.end(),
00525 hills_energy, hills_energy_gradients);
00526 new_hills_begin = hills.end();
00527
00528
00529
00530 if (comm == multiple_replicas) {
00531 for (size_t ir = 0; ir < replicas.size(); ir++) {
00532 replicas[ir]->project_hills (replicas[ir]->new_hills_begin,
00533 replicas[ir]->hills.end(),
00534 replicas[ir]->hills_energy,
00535 replicas[ir]->hills_energy_gradients);
00536 replicas[ir]->new_hills_begin = replicas[ir]->hills.end();
00537 }
00538 }
00539 }
00540
00541 std::vector<int> curr_bin = hills_energy->get_colvars_index();
00542 if (cvm::debug())
00543 cvm::log ("Metadynamics bias \""+this->name+"\""+
00544 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00545 ": current coordinates on the grid: "+
00546 cvm::to_str (curr_bin)+".\n");
00547
00548 if (hills_energy->index_ok (curr_bin)) {
00549
00550
00551
00552 bias_energy += hills_energy->value (curr_bin);
00553 for (size_t ic = 0; ic < colvars.size(); ic++) {
00554 cvm::real const *f = &(hills_energy_gradients->value (curr_bin));
00555 colvar_forces[ic].real_value += -1.0 * f[ic];
00556
00557 }
00558 if (comm == multiple_replicas)
00559 for (size_t ir = 0; ir < replicas.size(); ir++) {
00560 bias_energy += replicas[ir]->hills_energy->value (curr_bin);
00561 cvm::real const *f = &(replicas[ir]->hills_energy_gradients->value (curr_bin));
00562 for (size_t ic = 0; ic < colvars.size(); ic++) {
00563 colvar_forces[ic].real_value += -1.0 * f[ic];
00564 }
00565 }
00566
00567 } else {
00568
00569
00570
00571 calc_hills (hills_off_grid.begin(), hills_off_grid.end(), bias_energy);
00572 for (size_t ic = 0; ic < colvars.size(); ic++) {
00573 calc_hills_force (ic, hills_off_grid.begin(), hills_off_grid.end(), colvar_forces);
00574 }
00575
00576 if (comm == multiple_replicas)
00577 for (size_t ir = 0; ir < replicas.size(); ir++) {
00578 calc_hills (replicas[ir]->hills_off_grid.begin(),
00579 replicas[ir]->hills_off_grid.end(),
00580 bias_energy);
00581 for (size_t ic = 0; ic < colvars.size(); ic++) {
00582 calc_hills_force (ic,
00583 replicas[ir]->hills_off_grid.begin(),
00584 replicas[ir]->hills_off_grid.end(),
00585 colvar_forces);
00586 }
00587 }
00588 }
00589 }
00590
00591
00592
00593
00594 calc_hills (new_hills_begin, hills.end(), bias_energy);
00595 for (size_t ic = 0; ic < colvars.size(); ic++) {
00596 calc_hills_force (ic, new_hills_begin, hills.end(), colvar_forces);
00597 }
00598
00599 if (cvm::debug())
00600 cvm::log ("Hills energy = "+cvm::to_str (bias_energy)+
00601 ", hills forces = "+cvm::to_str (colvar_forces)+".\n");
00602
00603 if (cvm::debug())
00604 cvm::log ("Metadynamics bias \""+this->name+"\""+
00605 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00606 ": adding the forces from the other replicas.\n");
00607
00608 if (comm == multiple_replicas)
00609 for (size_t ir = 0; ir < replicas.size(); ir++) {
00610 calc_hills (replicas[ir]->new_hills_begin,
00611 replicas[ir]->hills.end(),
00612 bias_energy);
00613 for (size_t ic = 0; ic < colvars.size(); ic++) {
00614 calc_hills_force (ic,
00615 replicas[ir]->new_hills_begin,
00616 replicas[ir]->hills.end(),
00617 colvar_forces);
00618 }
00619 if (cvm::debug())
00620 cvm::log ("Hills energy = "+cvm::to_str (bias_energy)+
00621 ", hills forces = "+cvm::to_str (colvar_forces)+".\n");
00622 }
00623
00624 return bias_energy;
00625 }
00626
00627
00628 void colvarbias_meta::calc_hills (colvarbias_meta::hill_iter h_first,
00629 colvarbias_meta::hill_iter h_last,
00630 cvm::real &energy,
00631 std::vector<colvarvalue> const &colvar_values)
00632 {
00633 std::vector<colvarvalue> curr_values (colvars.size());
00634 for (size_t i = 0; i < colvars.size(); i++) {
00635 curr_values[i].type (colvars[i]->type());
00636 }
00637
00638 if (colvar_values.size()) {
00639 for (size_t i = 0; i < colvars.size(); i++) {
00640 curr_values[i] = colvar_values[i];
00641 }
00642 } else {
00643 for (size_t i = 0; i < colvars.size(); i++) {
00644 curr_values[i] = colvars[i]->value();
00645 }
00646 }
00647
00648 for (hill_iter h = h_first; h != h_last; h++) {
00649
00650
00651 cvm::real cv_sqdev = 0.0;
00652 for (size_t i = 0; i < colvars.size(); i++) {
00653 colvarvalue const &x = curr_values[i];
00654 colvarvalue const ¢er = h->centers[i];
00655 cvm::real const half_width = 0.5 * h->widths[i];
00656 cv_sqdev += (colvars[i]->dist2 (x, center)) / (half_width*half_width);
00657 }
00658
00659
00660 if (cv_sqdev > 23.0) {
00661
00662 h->value (0.0);
00663 } else {
00664 h->value (std::exp (-0.5*cv_sqdev));
00665 }
00666 energy += h->energy();
00667 }
00668 }
00669
00670
00671 void colvarbias_meta::calc_hills_force (size_t const &i,
00672 colvarbias_meta::hill_iter h_first,
00673 colvarbias_meta::hill_iter h_last,
00674 std::vector<colvarvalue> &forces,
00675 std::vector<colvarvalue> const &values)
00676 {
00677
00678 colvarvalue x (values.size() ? values[i].type() : colvars[i]->type());
00679 x = (values.size() ? values[i] : colvars[i]->value());
00680
00681
00682
00683
00684
00685 switch (colvars[i]->type()) {
00686
00687 case colvarvalue::type_scalar:
00688 for (hill_iter h = h_first; h != h_last; h++) {
00689 if (h->value() == 0.0) continue;
00690 colvarvalue const ¢er = h->centers[i];
00691 cvm::real const half_width = 0.5 * h->widths[i];
00692 forces[i].real_value +=
00693 ( h->weight() * h->value() * (0.5 / (half_width*half_width)) *
00694 (colvars[i]->dist2_lgrad (x, center)).real_value );
00695 }
00696 break;
00697
00698 case colvarvalue::type_vector:
00699 case colvarvalue::type_unitvector:
00700 for (hill_iter h = h_first; h != h_last; h++) {
00701 if (h->value() == 0.0) continue;
00702 colvarvalue const ¢er = h->centers[i];
00703 cvm::real const half_width = 0.5 * h->widths[i];
00704 forces[i].rvector_value +=
00705 ( h->weight() * h->value() * (0.5 / (half_width*half_width)) *
00706 (colvars[i]->dist2_lgrad (x, center)).rvector_value );
00707 }
00708 break;
00709
00710 case colvarvalue::type_quaternion:
00711 for (hill_iter h = h_first; h != h_last; h++) {
00712 if (h->value() == 0.0) continue;
00713 colvarvalue const ¢er = h->centers[i];
00714 cvm::real const half_width = 0.5 * h->widths[i];
00715 forces[i].quaternion_value +=
00716 ( h->weight() * h->value() * (0.5 / (half_width*half_width)) *
00717 (colvars[i]->dist2_lgrad (x, center)).quaternion_value );
00718 }
00719 break;
00720
00721 case colvarvalue::type_notset:
00722 break;
00723 }
00724 }
00725
00726
00727
00728
00729
00730
00731 void colvarbias_meta::project_hills (colvarbias_meta::hill_iter h_first,
00732 colvarbias_meta::hill_iter h_last,
00733 colvar_grid_scalar *he,
00734 colvar_grid_gradient *hg)
00735 {
00736 if (cvm::debug())
00737 cvm::log ("Metadynamics bias \""+this->name+"\""+
00738 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00739 ": projecting hills.\n");
00740
00741
00742
00743 std::vector<colvarvalue> colvar_values (colvars.size());
00744 std::vector<cvm::real> colvar_forces_scalar (colvars.size());
00745
00746 std::vector<int> he_ix = he->new_index();
00747 std::vector<int> hg_ix = (hg != NULL) ? hg->new_index() : std::vector<int> (0);
00748 cvm::real hills_energy_here = 0.0;
00749 std::vector<colvarvalue> hills_forces_here (colvars.size(), 0.0);
00750
00751 if (hg != NULL) {
00752
00753
00754 for ( ;
00755 (he->index_ok (he_ix)) && (hg->index_ok (hg_ix));
00756 ) {
00757
00758 for (size_t i = 0; i < colvars.size(); i++) {
00759 colvar_values[i] = hills_energy->bin_to_value_scalar (he_ix[i], i);
00760 }
00761
00762
00763 hills_energy_here = 0.0;
00764 calc_hills (h_first, h_last, hills_energy_here, colvar_values);
00765 he->acc_value (he_ix, hills_energy_here);
00766
00767 for (size_t i = 0; i < colvars.size(); i++) {
00768 hills_forces_here[i].reset();
00769 calc_hills_force (i, h_first, h_last, hills_forces_here, colvar_values);
00770 colvar_forces_scalar[i] = hills_forces_here[i].real_value;
00771 }
00772 hg->acc_force (hg_ix, &(colvar_forces_scalar.front()));
00773
00774 he->incr (he_ix);
00775 hg->incr (hg_ix);
00776 }
00777
00778 } else {
00779
00780
00781
00782 for ( ; (he->index_ok (he_ix)); ) {
00783
00784 for (size_t i = 0; i < colvars.size(); i++) {
00785 colvar_values[i] = hills_energy->bin_to_value_scalar (he_ix[i], i);
00786 }
00787
00788 hills_energy_here = 0.0;
00789 calc_hills (h_first, h_last, hills_energy_here, colvar_values);
00790 he->acc_value (he_ix, hills_energy_here);
00791
00792 he->incr (he_ix);
00793 }
00794 }
00795
00796 if (! keep_hills) {
00797 hills.erase (hills.begin(), hills.end());
00798 }
00799 }
00800
00801
00802 void colvarbias_meta::recount_hills_off_grid (colvarbias_meta::hill_iter h_first,
00803 colvarbias_meta::hill_iter h_last,
00804 colvar_grid_scalar *he)
00805 {
00806 hills_off_grid.clear();
00807
00808 for (hill_iter h = h_first; h != h_last; h++) {
00809 cvm::real const min_dist = hills_energy->bin_distance_from_boundaries (h->centers);
00810 if (min_dist < (3.0 * std::floor (hill_width)) + 1.0) {
00811 hills_off_grid.push_back (*h);
00812 }
00813 }
00814 }
00815
00816
00817
00818
00819
00820
00821
00822
00823 void colvarbias_meta::update_replicas_registry()
00824 {
00825 if (cvm::debug())
00826 cvm::log ("Metadynamics bias \""+this->name+"\""+
00827 ": updating the list of replicas, currently containing "+
00828 cvm::to_str (replicas.size())+" elements.\n");
00829
00830 {
00831
00832 std::string line ("");
00833 std::ifstream reg_file (replicas_registry_file.c_str());
00834 if (reg_file.good()) {
00835 replicas_registry.clear();
00836 while (colvarparse::getline_nocomments (reg_file, line))
00837 replicas_registry.append (line+"\n");
00838 } else {
00839 cvm::fatal_error ("Error: failed to open file \""+replicas_registry_file+
00840 "\" for reading.\n");
00841 }
00842 }
00843
00844
00845 std::istringstream reg_is (replicas_registry);
00846 if (reg_is.good()) {
00847
00848 std::string new_replica ("");
00849 std::string new_replica_file ("");
00850 while ((reg_is >> new_replica) && new_replica.size() &&
00851 (reg_is >> new_replica_file) && new_replica_file.size()) {
00852
00853 if (new_replica == this->replica_id) {
00854
00855 if (cvm::debug())
00856 cvm::log ("Metadynamics bias \""+this->name+"\""+
00857 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00858 ": skipping this replica's own record: \""+
00859 new_replica+"\", \""+new_replica_file+"\"\n");
00860 new_replica_file.clear();
00861 new_replica.clear();
00862 continue;
00863 }
00864
00865 bool already_loaded = false;
00866 for (size_t ir = 0; ir < replicas.size(); ir++) {
00867 if (new_replica == (replicas[ir])->replica_id) {
00868
00869 if (cvm::debug())
00870 cvm::log ("Metadynamics bias \""+this->name+"\""+
00871 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
00872 ": skipping a replica already loaded, \""+
00873 (replicas[ir])->replica_id+"\".\n");
00874 already_loaded = true;
00875 break;
00876 }
00877 }
00878
00879 if (!already_loaded) {
00880
00881 cvm::log ("Metadynamics bias \""+this->name+"\""+
00882 ": accessing replica \""+new_replica+"\".\n");
00883 replicas.push_back (new colvarbias_meta());
00884 (replicas.back())->replica_id = new_replica;
00885 (replicas.back())->replica_list_file = new_replica_file;
00886 (replicas.back())->replica_state_file = "";
00887 (replicas.back())->replica_state_file_in_sync = false;
00888
00889
00890 (replicas.back())->name = this->name;
00891 (replicas.back())->colvars = colvars;
00892 (replicas.back())->use_grids = use_grids;
00893 (replicas.back())->dump_fes = false;
00894 (replicas.back())->expand_grids = false;
00895 (replicas.back())->rebin_grids = false;
00896 (replicas.back())->keep_hills = false;
00897 (replicas.back())->colvar_forces = colvar_forces;
00898
00899 (replicas.back())->comm = multiple_replicas;
00900
00901 if (use_grids) {
00902 (replicas.back())->hills_energy = new colvar_grid_scalar (colvars);
00903 (replicas.back())->hills_energy_gradients = new colvar_grid_gradient (colvars);
00904 }
00905 }
00906 }
00907
00908
00909 new_replica_file = "";
00910 new_replica = "";
00911 } else {
00912 cvm::fatal_error ("Error: cannot read the replicas registry file \""+
00913 replicas_registry+"\".\n");
00914 }
00915
00916
00917 for (size_t ir = 0; ir < replicas.size(); ir++) {
00918 if (cvm::debug())
00919 cvm::log ("Metadynamics bias \""+this->name+"\""+
00920 ": reading the list file for replica \""+
00921 (replicas[ir])->replica_id+"\".\n");
00922
00923 std::ifstream list_is ((replicas[ir])->replica_list_file.c_str());
00924 std::string key;
00925 std::string new_state_file, new_hills_file;
00926 if (!(list_is >> key) ||
00927 !(key == std::string ("stateFile")) ||
00928 !(list_is >> new_state_file) ||
00929 !(list_is >> key) ||
00930 !(key == std::string ("hillsFile")) ||
00931 !(list_is >> new_hills_file)) {
00932 cvm::log ("Metadynamics bias \""+this->name+"\""+
00933 ": failed to read the file \""+
00934 (replicas[ir])->replica_list_file+"\": will try again after "+
00935 cvm::to_str (replica_update_freq)+" steps.\n");
00936 (replicas[ir])->update_status++;
00937 } else {
00938 (replicas[ir])->update_status = 0;
00939 if (new_state_file != (replicas[ir])->replica_state_file) {
00940 cvm::log ("Metadynamics bias \""+this->name+"\""+
00941 ": replica \""+(replicas[ir])->replica_id+
00942 "\" has supplied a new state file, \""+new_state_file+
00943 "\".\n");
00944 (replicas[ir])->replica_state_file_in_sync = false;
00945 (replicas[ir])->replica_state_file = new_state_file;
00946 (replicas[ir])->replica_hills_file = new_hills_file;
00947 }
00948 }
00949 }
00950
00951
00952 if (cvm::debug())
00953 cvm::log ("Metadynamics bias \""+this->name+"\": the list of replicas contains "+
00954 cvm::to_str (replicas.size())+" elements.\n");
00955 }
00956
00957
00958 void colvarbias_meta::read_replica_files()
00959 {
00960 for (size_t ir = 0; ir < replicas.size(); ir++) {
00961
00962 if (! (replicas[ir])->replica_state_file_in_sync) {
00963
00964 (replicas[ir])->replica_hills_file_pos = 0;
00965 }
00966
00967
00968 if ( (! (replicas[ir])->has_data) ||
00969 (! (replicas[ir])->replica_state_file_in_sync) ) {
00970
00971 cvm::log ("Metadynamics bias \""+this->name+"\""+
00972 ": reading the state of replica \""+
00973 (replicas[ir])->replica_id+"\" from file \""+
00974 (replicas[ir])->replica_state_file+"\".\n");
00975
00976 std::ifstream is ((replicas[ir])->replica_state_file.c_str());
00977 if (! (replicas[ir])->read_restart (is)) {
00978 cvm::log ("Reading from file \""+(replicas[ir])->replica_state_file+
00979 "\" failed or incomplete: will try again in "+
00980 cvm::to_str (replica_update_freq)+" steps.\n");
00981 } else {
00982
00983 (replicas[ir])->replica_state_file_in_sync = true;
00984 (replicas[ir])->update_status = 0;
00985 }
00986 is.close();
00987 }
00988
00989
00990 if ((replicas[ir])->replica_hills_file.size()) {
00991
00992 if (cvm::debug())
00993 cvm::log ("Metadynamics bias \""+this->name+"\""+
00994 ": checking for new hills from replica \""+
00995 (replicas[ir])->replica_id+"\" in the file \""+
00996 (replicas[ir])->replica_hills_file+"\".\n");
00997
00998
00999
01000
01001 std::ifstream is ((replicas[ir])->replica_hills_file.c_str());
01002 if (is.good()) {
01003
01004
01005 is.seekg ((replicas[ir])->replica_hills_file_pos, std::ios::beg);
01006 if (!is.good()){
01007
01008
01009 is.clear();
01010 is.seekg (0, std::ios::beg);
01011
01012 (replicas[ir])->replica_hills_file_pos = 0;
01013
01014 (replicas[ir])->replica_state_file_in_sync = false;
01015
01016 (replicas[ir])->update_status++;
01017 cvm::log ("Failed to read the file \""+(replicas[ir])->replica_hills_file+
01018 "\" at the previous position: will try again in "+
01019 cvm::to_str (replica_update_freq)+" steps.\n");
01020 } else {
01021
01022 while ((replicas[ir])->read_hill (is)) {
01023
01024 cvm::log ("Metadynamics bias \""+this->name+"\""+
01025 ": received a hill from replica \""+
01026 (replicas[ir])->replica_id+
01027 "\" at step "+
01028 cvm::to_str (((replicas[ir])->hills.back()).it)+".\n");
01029 }
01030 is.clear();
01031
01032 (replicas[ir])->replica_hills_file_pos = is.tellg();
01033 if (cvm::debug())
01034 cvm::log ("Metadynamics bias \""+this->name+"\""+
01035 ": stopped reading file \""+(replicas[ir])->replica_hills_file+
01036 "\" at position "+
01037 cvm::to_str ((replicas[ir])->replica_hills_file_pos)+".\n");
01038
01039
01040 is.seekg (0, std::ios::end);
01041 if (is.tellg() > (replicas[ir])->replica_hills_file_pos+1) {
01042 (replicas[ir])->update_status++;
01043 } else {
01044 (replicas[ir])->update_status = 0;
01045 }
01046 }
01047
01048 } else {
01049 cvm::fatal_error ("Error: cannot read from file \""+
01050 (replicas[ir])->replica_hills_file+"\".\n");
01051 }
01052 is.close();
01053 }
01054
01055 size_t const n_flush = (replica_update_freq/new_hill_freq + 1);
01056 if ((replicas[ir])->update_status > 3*n_flush) {
01057
01058 cvm::log ("Warning: in metadynamics bias \""+this->name+"\""+
01059 ": failet do read completely output files from replica \""+
01060 (replicas[ir])->replica_id+
01061 "\" after more than "+
01062 cvm::to_str (n_flush*new_hill_freq)+
01063 " steps. Please check that that simulation is still running.\n");
01064 }
01065 }
01066 }
01067
01068
01069
01070
01071
01072
01073
01074 std::istream & colvarbias_meta::read_restart (std::istream& is)
01075 {
01076 size_t const start_pos = is.tellg();
01077
01078 if (comm == single_replica) {
01079
01080
01081 cvm::log ("Restarting metadynamics bias \""+this->name+"\""+
01082 ".\n");
01083 }
01084 std::string key, brace, conf;
01085 if ( !(is >> key) || !(key == "metadynamics") ||
01086 !(is >> brace) || !(brace == "{") ||
01087 !(is >> colvarparse::read_block ("configuration", conf)) ) {
01088
01089 if (comm == single_replica)
01090 cvm::log ("Error: in reading restart configuration for metadynamics bias \""+
01091 this->name+"\""+
01092 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
01093 (replica_state_file_in_sync ? ("at position "+
01094 cvm::to_str (start_pos)+
01095 " in the state file") : "")+".\n");
01096 is.clear();
01097 is.seekg (start_pos, std::ios::beg);
01098 is.setstate (std::ios::failbit);
01099 return is;
01100 }
01101
01102 std::string name = "";
01103 if ( colvarparse::get_keyval (conf, "name", name,
01104 std::string (""), colvarparse::parse_silent) &&
01105 (name != this->name) )
01106 cvm::fatal_error ("Error: in the restart file, the "
01107 "\"metadynamics\" block has a different name: different system?\n");
01108
01109 if (name.size() == 0) {
01110 cvm::fatal_error ("Error: \"metadynamics\" block within the restart file "
01111 "has no identifiers.\n");
01112 }
01113
01114 if (comm != single_replica) {
01115 std::string replica = "";
01116 if (colvarparse::get_keyval (conf, "replicaID", replica,
01117 std::string (""), colvarparse::parse_silent) &&
01118 (replica != this->replica_id))
01119 cvm::fatal_error ("Error: in the restart file, the "
01120 "\"metadynamics\" block has a different replicaID: different system?\n");
01121
01122 colvarparse::get_keyval (conf, "step", state_file_step,
01123 cvm::step_absolute(), colvarparse::parse_silent);
01124 }
01125
01126 bool grids_from_restart_file = use_grids;
01127
01128 if (use_grids) {
01129
01130 if (expand_grids) {
01131
01132
01133
01134
01135 delete hills_energy;
01136 delete hills_energy_gradients;
01137 hills_energy = new colvar_grid_scalar (colvars);
01138 hills_energy_gradients = new colvar_grid_gradient (colvars);
01139 }
01140
01141 colvar_grid_scalar *hills_energy_backup = NULL;
01142 colvar_grid_gradient *hills_energy_gradients_backup = NULL;
01143
01144 if (has_data) {
01145 if (cvm::debug())
01146 cvm::log ("Backupping grids for metadynamics bias \""+
01147 this->name+"\""+
01148 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+".\n");
01149 hills_energy_backup = hills_energy;
01150 hills_energy_gradients_backup = hills_energy_gradients;
01151 hills_energy = new colvar_grid_scalar (colvars);
01152 hills_energy_gradients = new colvar_grid_gradient (colvars);
01153 }
01154
01155 size_t const hills_energy_pos = is.tellg();
01156 if (!(is >> key)) {
01157 if (hills_energy_backup != NULL) {
01158 delete hills_energy;
01159 delete hills_energy_gradients;
01160 hills_energy = hills_energy_backup;
01161 hills_energy_gradients = hills_energy_gradients_backup;
01162 }
01163 is.clear();
01164 is.seekg (hills_energy_pos, std::ios::beg);
01165 is.setstate (std::ios::failbit);
01166 return is;
01167 } else if (!(key == std::string ("hills_energy")) ||
01168 !(hills_energy->read_restart (is))) {
01169 is.clear();
01170 is.seekg (hills_energy_pos, std::ios::beg);
01171 grids_from_restart_file = false;
01172 if (!rebin_grids) {
01173 if (hills_energy_backup == NULL)
01174 cvm::fatal_error ("Error: couldn't read the free energy grid for metadynamics bias \""+
01175 this->name+"\""+
01176 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
01177 "; if useGrids was off when the state file was written, "
01178 "enable rebinGrids now to regenerate the grids.\n");
01179 else {
01180 if (comm == single_replica)
01181 cvm::log ("Error: couldn't read the free energy grid for metadynamics bias \""+
01182 this->name+"\".\n");
01183 delete hills_energy;
01184 delete hills_energy_gradients;
01185 hills_energy = hills_energy_backup;
01186 hills_energy_gradients = hills_energy_gradients_backup;
01187 is.setstate (std::ios::failbit);
01188 return is;
01189 }
01190 }
01191 }
01192
01193 size_t const hills_energy_gradients_pos = is.tellg();
01194 if (!(is >> key)) {
01195 if (hills_energy_backup != NULL) {
01196 delete hills_energy;
01197 delete hills_energy_gradients;
01198 hills_energy = hills_energy_backup;
01199 hills_energy_gradients = hills_energy_gradients_backup;
01200 }
01201 is.clear();
01202 is.seekg (hills_energy_gradients_pos, std::ios::beg);
01203 is.setstate (std::ios::failbit);
01204 return is;
01205 } else if (!(key == std::string ("hills_energy_gradients")) ||
01206 !(hills_energy_gradients->read_restart (is))) {
01207 is.clear();
01208 is.seekg (hills_energy_gradients_pos, std::ios::beg);
01209 grids_from_restart_file = false;
01210 if (hills_energy_backup == NULL) {
01211 if (!rebin_grids)
01212 cvm::fatal_error ("Error: couldn't read the free energy gradients grid for metadynamics bias \""+
01213 this->name+"\""+
01214 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
01215 "; if useGrids was off when the state file was written, "
01216 "enable rebinGrids now to regenerate the grids.\n");
01217 else {
01218 if (comm == single_replica)
01219 cvm::log ("Error: couldn't read the free energy gradients grid for metadynamics bias \""+
01220 this->name+"\".\n");
01221 delete hills_energy;
01222 delete hills_energy_gradients;
01223 hills_energy = hills_energy_backup;
01224 hills_energy_gradients = hills_energy_gradients_backup;
01225 is.setstate (std::ios::failbit);
01226 return is;
01227 }
01228 }
01229 }
01230
01231 if (cvm::debug())
01232 cvm::log ("Successfully read new grids for bias \""+
01233 this->name+"\""+
01234 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+"\n");
01235
01236 if (hills_energy_backup != NULL) {
01237
01238
01239 if (cvm::debug())
01240 cvm::log ("Deallocating the older grids.\n");
01241
01242 delete hills_energy_backup;
01243 delete hills_energy_gradients_backup;
01244 }
01245 }
01246
01247 bool const existing_hills = (hills.size() > 0);
01248 size_t const old_hills_size = hills.size();
01249 hill_iter old_hills_end = hills.end();
01250 hill_iter old_hills_off_grid_end = hills_off_grid.end();
01251
01252
01253 while (read_hill (is)) {
01254 if (cvm::debug())
01255 cvm::log ("Read a previously saved hill under the "
01256 "metadynamics bias \""+
01257 this->name+"\", created at step "+
01258 cvm::to_str ((hills.back()).it)+".\n");
01259 }
01260 is.clear();
01261 new_hills_begin = hills.end();
01262 if (grids_from_restart_file) {
01263 if (hills.size() > old_hills_size)
01264 cvm::log ("Read "+cvm::to_str (hills.size())+
01265 " hills in addition to the grids.\n");
01266 } else {
01267 if (hills.size())
01268 cvm::log ("Read "+cvm::to_str (hills.size())+" hills.\n");
01269 }
01270
01271 if (rebin_grids) {
01272
01273
01274
01275
01276
01277 colvar_grid_scalar *new_hills_energy =
01278 new colvar_grid_scalar (colvars);
01279 colvar_grid_gradient *new_hills_energy_gradients =
01280 new colvar_grid_gradient (colvars);
01281
01282 if (!grids_from_restart_file || (keep_hills && hills.size())) {
01283
01284 cvm::log ("Rebinning the energy and forces grids from "+
01285 cvm::to_str (hills.size())+" hills (this may take a while)...\n");
01286 project_hills (hills.begin(), hills.end(),
01287 new_hills_energy, new_hills_energy_gradients);
01288 cvm::log ("rebinning done.\n");
01289
01290 } else {
01291
01292 cvm::log ("Rebinning the energy and forces grids "
01293 "from the grids in the restart file.\n");
01294 new_hills_energy->map_grid (*hills_energy);
01295 new_hills_energy_gradients->map_grid (*hills_energy_gradients);
01296 }
01297
01298 delete hills_energy;
01299 delete hills_energy_gradients;
01300 hills_energy = new_hills_energy;
01301 hills_energy_gradients = new_hills_energy_gradients;
01302
01303
01304
01305 if (hills.size())
01306 recount_hills_off_grid (hills.begin(), hills.end(), hills_energy);
01307 }
01308
01309 if (use_grids) {
01310 if (hills_off_grid.size()) {
01311 cvm::log (cvm::to_str (hills_off_grid.size())+" hills are near the "
01312 "grid boundaries: they will be computed analytically "
01313 "and saved to the state files.\n");
01314 }
01315 }
01316
01317 is >> brace;
01318 if (brace != "}") {
01319 cvm::log ("Incomplete restart information for metadynamics bias \""+
01320 this->name+"\""+
01321 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+
01322 ": no closing brace at position "+
01323 cvm::to_str (is.tellg())+" in the file.\n");
01324 is.setstate (std::ios::failbit);
01325 return is;
01326 }
01327
01328 if (cvm::debug())
01329 cvm::log ("colvarbias_meta::read_restart() done\n");
01330
01331 if (existing_hills) {
01332 hills.erase (hills.begin(), old_hills_end);
01333 hills_off_grid.erase (hills_off_grid.begin(), old_hills_off_grid_end);
01334 }
01335
01336 has_data = true;
01337
01338 if (comm != single_replica) {
01339 read_replica_files();
01340 }
01341
01342 return is;
01343 }
01344
01345
01346 std::istream & colvarbias_meta::read_hill (std::istream &is)
01347 {
01348 if (!is) return is;
01349
01350 size_t const start_pos = is.tellg();
01351
01352 std::string data;
01353 if ( !(is >> read_block ("hill", data)) ) {
01354 is.clear();
01355 is.seekg (start_pos, std::ios::beg);
01356 is.setstate (std::ios::failbit);
01357 return is;
01358 }
01359
01360 size_t h_it;
01361 get_keyval (data, "step", h_it, 0, parse_silent);
01362 if (h_it <= state_file_step) {
01363 if (cvm::debug())
01364 cvm::log ("Skipping a hill older than the state file for metadynamics bias \""+
01365 this->name+"\""+
01366 ((comm != single_replica) ? ", replica \""+replica_id+"\"" : "")+"\n");
01367 return is;
01368 }
01369
01370 cvm::real h_weight;
01371 get_keyval (data, "weight", h_weight, hill_weight, parse_silent);
01372
01373 std::vector<colvarvalue> h_centers (colvars.size());
01374 for (size_t i = 0; i < colvars.size(); i++) {
01375 h_centers[i].type ((colvars[i]->value()).type());
01376 }
01377 {
01378
01379
01380 std::string centers_input;
01381 key_lookup (data, "centers", centers_input);
01382 std::istringstream centers_is (centers_input);
01383 for (size_t i = 0; i < colvars.size(); i++) {
01384 centers_is >> h_centers[i];
01385 }
01386 }
01387
01388 std::vector<cvm::real> h_widths (colvars.size());
01389 get_keyval (data, "widths", h_widths,
01390 std::vector<cvm::real> (colvars.size(), (std::sqrt (2.0 * PI) / 2.0)),
01391 parse_silent);
01392
01393 std::string h_replica = "";
01394 if (comm != single_replica) {
01395 get_keyval (data, "replicaID", h_replica, replica_id, parse_silent);
01396 if (h_replica != replica_id)
01397 cvm::fatal_error ("Error: trying to read a hill created by replica \""+h_replica+
01398 "\" for replica \""+replica_id+
01399 "\"; did you swap output files?\n");
01400 }
01401
01402 hill_iter const hills_end = hills.end();
01403 hills.push_back (hill (h_it, h_weight, h_centers, h_widths, h_replica));
01404 if (new_hills_begin == hills_end) {
01405
01406 new_hills_begin = hills.end();
01407 new_hills_begin--;
01408 }
01409
01410 if (use_grids) {
01411
01412
01413 cvm::real const min_dist =
01414 hills_energy->bin_distance_from_boundaries ((hills.back()).centers);
01415 if (min_dist < (3.0 * std::floor (hill_width)) + 1.0) {
01416 hills_off_grid.push_back (hills.back());
01417 }
01418 }
01419
01420 has_data = true;
01421 return is;
01422 }
01423
01424
01425
01426
01427
01428
01429
01430
01431 std::ostream & colvarbias_meta::write_restart (std::ostream& os)
01432 {
01433 os << "metadynamics {\n"
01434 << " configuration {\n"
01435 << " step " << cvm::step_absolute() << "\n"
01436 << " name " << this->name << "\n";
01437 if (this->comm != single_replica)
01438 os << " replicaID " << this->replica_id << "\n";
01439 os << " }\n\n";
01440
01441 if (use_grids) {
01442
01443
01444
01445 project_hills (new_hills_begin, hills.end(),
01446 hills_energy, hills_energy_gradients);
01447 new_hills_begin = hills.end();
01448
01449
01450 os << " hills_energy\n";
01451 hills_energy->write_restart (os);
01452 os << " hills_energy_gradients\n";
01453 hills_energy_gradients->write_restart (os);
01454 }
01455
01456 if ( (!use_grids) || keep_hills ) {
01457
01458 for (std::list<hill>::const_iterator h = this->hills.begin();
01459 h != this->hills.end();
01460 h++) {
01461 os << *h;
01462 }
01463 } else {
01464
01465 for (std::list<hill>::const_iterator h = this->hills_off_grid.begin();
01466 h != this->hills_off_grid.end();
01467 h++) {
01468 os << *h;
01469 }
01470 }
01471
01472 os << "}\n\n";
01473
01474 if (comm != single_replica) {
01475 write_replica_state_file();
01476
01477
01478 for (size_t ir = 0; ir < replicas.size(); ir++) {
01479 (replicas[ir])->replica_state_file_in_sync = false;
01480 }
01481 }
01482
01483 if (dump_fes) {
01484 write_pmf();
01485 }
01486
01487 return os;
01488 }
01489
01490
01491 void colvarbias_meta::write_pmf()
01492 {
01493
01494 colvar_grid_scalar *pmf = new colvar_grid_scalar (colvars);
01495
01496 std::string fes_file_name_prefix (cvm::output_prefix);
01497
01498 if ((cvm::n_meta_biases > 1) || (cvm::n_abf_biases > 0)) {
01499
01500
01501
01502
01503 fes_file_name_prefix += ("."+this->name);
01504 }
01505
01506 if ((comm == single_replica) || (dump_replica_fes)) {
01507
01508 pmf->reset();
01509 pmf->add_grid (*hills_energy);
01510 cvm::real const max = pmf->maximum_value();
01511 pmf->add_constant (-1.0 * max);
01512 if (!well_tempered) {
01513 pmf->multiply_constant (-1.0);
01514 } else {
01515 cvm::real well_temper_scale = (bias_temperature + cvm::temperature()) / bias_temperature;
01516 pmf->multiply_constant (-well_temper_scale);
01517 }
01518
01519 std::string const fes_file_name (fes_file_name_prefix +
01520 ((comm != single_replica) ? ".partial" : "") +
01521 (dump_fes_save ?
01522 "."+cvm::to_str (cvm::step_absolute()) : "") +
01523 ".pmf");
01524 cvm::backup_file (fes_file_name.c_str());
01525 std::ofstream fes_os (fes_file_name.c_str());
01526 pmf->write_multicol (fes_os);
01527 fes_os.close();
01528 }
01529
01530 if (comm != single_replica) {
01531
01532 pmf->reset();
01533 pmf->add_grid (*hills_energy);
01534 for (size_t ir = 0; ir < replicas.size(); ir++) {
01535 pmf->add_grid (*(replicas[ir]->hills_energy));
01536 }
01537 cvm::real const max = pmf->maximum_value();
01538 pmf->add_constant (-1.0 * max);
01539 pmf->multiply_constant (-1.0);
01540 std::string const fes_file_name (fes_file_name_prefix +
01541 (dump_fes_save ?
01542 "."+cvm::to_str (cvm::step_absolute()) : "") +
01543 ".pmf");
01544 cvm::backup_file (fes_file_name.c_str());
01545 std::ofstream fes_os (fes_file_name.c_str());
01546 pmf->write_multicol (fes_os);
01547 fes_os.close();
01548 }
01549
01550 delete pmf;
01551 }
01552
01553
01554
01555 void colvarbias_meta::write_replica_state_file()
01556 {
01557
01558
01559 cvm::backup_file (replica_state_file.c_str());
01560 std::ofstream rep_state_os (replica_state_file.c_str());
01561 if (!rep_state_os.good())
01562 cvm::fatal_error ("Error: in opening file \""+
01563 replica_state_file+"\" for writing.\n");
01564
01565 rep_state_os.setf (std::ios::scientific, std::ios::floatfield);
01566 rep_state_os << "\n"
01567 << "metadynamics {\n"
01568 << " configuration {\n"
01569 << " name " << this->name << "\n"
01570 << " step " << cvm::step_absolute() << "\n";
01571 if (this->comm != single_replica) {
01572 rep_state_os << " replicaID " << this->replica_id << "\n";
01573 }
01574 rep_state_os << " }\n\n";
01575 rep_state_os << " hills_energy\n";
01576 rep_state_os << std::setprecision (cvm::cv_prec)
01577 << std::setw (cvm::cv_width);
01578 hills_energy->write_restart (rep_state_os);
01579 rep_state_os << " hills_energy_gradients\n";
01580 rep_state_os << std::setprecision (cvm::cv_prec)
01581 << std::setw (cvm::cv_width);
01582 hills_energy_gradients->write_restart (rep_state_os);
01583
01584 if ( (!use_grids) || keep_hills ) {
01585
01586 for (std::list<hill>::const_iterator h = this->hills.begin();
01587 h != this->hills.end();
01588 h++) {
01589 rep_state_os << *h;
01590 }
01591 } else {
01592
01593 for (std::list<hill>::const_iterator h = this->hills_off_grid.begin();
01594 h != this->hills_off_grid.end();
01595 h++) {
01596 rep_state_os << *h;
01597 }
01598 }
01599
01600 rep_state_os << "}\n\n";
01601 rep_state_os.close();
01602
01603
01604 replica_hills_os.close();
01605 replica_hills_os.open (replica_hills_file.c_str());
01606 if (!replica_hills_os.good())
01607 cvm::fatal_error ("Error: in opening file \""+
01608 replica_hills_file+"\" for writing.\n");
01609 replica_hills_os.setf (std::ios::scientific, std::ios::floatfield);
01610 }
01611
01612 std::string colvarbias_meta::hill::output_traj()
01613 {
01614 std::ostringstream os;
01615 os.setf (std::ios::fixed, std::ios::floatfield);
01616 os << std::setw (cvm::it_width) << it << " ";
01617
01618 os.setf (std::ios::scientific, std::ios::floatfield);
01619
01620 os << " ";
01621 for (size_t i = 0; i < centers.size(); i++) {
01622 os << " ";
01623 os << std::setprecision (cvm::cv_prec)
01624 << std::setw (cvm::cv_width) << centers[i];
01625 }
01626
01627 os << " ";
01628 for (size_t i = 0; i < widths.size(); i++) {
01629 os << " ";
01630 os << std::setprecision (cvm::cv_prec)
01631 << std::setw (cvm::cv_width) << widths[i];
01632 }
01633
01634 os << " ";
01635 os << std::setprecision (cvm::en_prec)
01636 << std::setw (cvm::en_width) << W << "\n";
01637
01638 return os.str();
01639 }
01640
01641
01642 std::ostream & operator << (std::ostream &os, colvarbias_meta::hill const &h)
01643 {
01644 os.setf (std::ios::scientific, std::ios::floatfield);
01645
01646 os << "hill {\n";
01647 os << " step " << std::setw (cvm::it_width) << h.it << "\n";
01648 os << " weight "
01649 << std::setprecision (cvm::en_prec)
01650 << std::setw (cvm::en_width)
01651 << h.W << "\n";
01652
01653 if (h.replica.size())
01654 os << " replicaID " << h.replica << "\n";
01655
01656 os << " centers ";
01657 for (size_t i = 0; i < (h.centers).size(); i++) {
01658 os << " "
01659 << std::setprecision (cvm::cv_prec)
01660 << std::setw (cvm::cv_width)
01661 << h.centers[i];
01662 }
01663 os << "\n";
01664
01665 os << " widths ";
01666 for (size_t i = 0; i < (h.widths).size(); i++) {
01667 os << " "
01668 << std::setprecision (cvm::cv_prec)
01669 << std::setw (cvm::cv_width)
01670 << h.widths[i];
01671 }
01672 os << "\n";
01673
01674 os << "}\n";
01675
01676 return os;
01677 }