1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670
|
// -*- c++ -*-
// This file is part of the Collective Variables module (Colvars).
// The original version of Colvars and its updates are located at:
// https://github.com/colvars/colvars
// Please update all Colvars source files before making any changes.
// If you wish to distribute your changes, please submit them to the
// Colvars repository at GitHub.
#include "colvarmodule.h"
#include "colvarproxy.h"
#include "colvarvalue.h"
#include "colvarbias.h"
#include "colvargrid.h"
colvarbias::colvarbias(char const *key)
: bias_type(to_lower_cppstr(key))
{
init_cvb_requires();
rank = 1;
has_data = false;
b_output_energy = false;
reset();
state_file_step = 0;
description = "uninitialized " + cvm::to_str(key) + " bias";
}
int colvarbias::init(std::string const &conf)
{
colvarparse::init(conf);
size_t i = 0;
if (name.size() == 0) {
// first initialization
cvm::log("Initializing a new \""+bias_type+"\" instance.\n");
rank = cvm::main()->num_biases_type(bias_type);
get_keyval(conf, "name", name, bias_type+cvm::to_str(rank));
{
colvarbias *bias_with_name = cvm::bias_by_name(this->name);
if (bias_with_name != NULL) {
if ((bias_with_name->rank != this->rank) ||
(bias_with_name->bias_type != this->bias_type)) {
cvm::error("Error: this bias cannot have the same name, \""+this->name+
"\", as another bias.\n", INPUT_ERROR);
return INPUT_ERROR;
}
}
}
description = "bias " + name;
{
// lookup the associated colvars
std::vector<std::string> colvar_names;
if (get_keyval(conf, "colvars", colvar_names)) {
if (num_variables()) {
cvm::error("Error: cannot redefine the colvars that a bias was already defined on.\n",
INPUT_ERROR);
return INPUT_ERROR;
}
for (i = 0; i < colvar_names.size(); i++) {
add_colvar(colvar_names[i]);
}
}
}
if (!num_variables()) {
cvm::error("Error: no collective variables specified.\n", INPUT_ERROR);
return INPUT_ERROR;
}
} else {
cvm::log("Reinitializing bias \""+name+"\".\n");
}
output_prefix = cvm::output_prefix();
get_keyval(conf, "outputEnergy", b_output_energy, b_output_energy);
get_keyval(conf, "timeStepFactor", time_step_factor, 1);
if (time_step_factor < 1) {
cvm::error("Error: timeStepFactor must be 1 or greater.\n");
return COLVARS_ERROR;
}
// Now that children are defined, we can solve dependencies
enable(f_cvb_active);
if (cvm::debug()) print_state();
return COLVARS_OK;
}
int colvarbias::reset()
{
bias_energy = 0.0;
for (size_t i = 0; i < num_variables(); i++) {
colvar_forces[i].reset();
}
return COLVARS_OK;
}
colvarbias::colvarbias()
: colvarparse(), has_data(false)
{}
colvarbias::~colvarbias()
{
colvarbias::clear();
}
int colvarbias::clear()
{
free_children_deps();
// Remove references to this bias from colvars
for (std::vector<colvar *>::iterator cvi = colvars.begin();
cvi != colvars.end();
++cvi) {
for (std::vector<colvarbias *>::iterator bi = (*cvi)->biases.begin();
bi != (*cvi)->biases.end();
++bi) {
if ( *bi == this) {
(*cvi)->biases.erase(bi);
break;
}
}
}
colvarmodule *cv = cvm::main();
// ...and from the colvars module
for (std::vector<colvarbias *>::iterator bi = cv->biases.begin();
bi != cv->biases.end();
++bi) {
if ( *bi == this) {
cv->biases.erase(bi);
break;
}
}
return COLVARS_OK;
}
int colvarbias::clear_state_data()
{
// no mutable content to delete for base class
return COLVARS_OK;
}
int colvarbias::add_colvar(std::string const &cv_name)
{
if (colvar *cv = cvm::colvar_by_name(cv_name)) {
if (cvm::debug()) {
cvm::log("Applying this bias to collective variable \""+
cv->name+"\".\n");
}
colvars.push_back(cv);
colvar_forces.push_back(colvarvalue());
colvar_forces.back().type(cv->value()); // make sure each force is initialized to zero
colvar_forces.back().is_derivative(); // colvar constraints are not applied to the force
colvar_forces.back().reset();
previous_colvar_forces.push_back(colvar_forces.back());
cv->biases.push_back(this); // add back-reference to this bias to colvar
if (is_enabled(f_cvb_apply_force)) {
cv->enable(f_cv_gradient);
}
// Add dependency link.
// All biases need at least the value of each colvar
// although possibly not at all timesteps
add_child(cv);
} else {
cvm::error("Error: cannot find a colvar named \""+
cv_name+"\".\n", INPUT_ERROR);
return INPUT_ERROR;
}
return COLVARS_OK;
}
int colvarbias::update()
{
if (cvm::debug()) {
cvm::log("Updating the "+bias_type+" bias \""+this->name+"\".\n");
}
has_data = true;
bias_energy = 0.0;
for (size_t ir = 0; ir < num_variables(); ir++) {
colvar_forces[ir].reset();
}
return COLVARS_OK;
}
void colvarbias::communicate_forces()
{
size_t i = 0;
for (i = 0; i < num_variables(); i++) {
if (cvm::debug()) {
cvm::log("Communicating a force to colvar \""+
variables(i)->name+"\".\n");
}
// Impulse-style multiple timestep
// Note that biases with different values of time_step_factor
// may send forces to the same colvar
// which is why rescaling has to happen now: the colvar is not
// aware of this bias' time_step_factor
variables(i)->add_bias_force(cvm::real(time_step_factor) * colvar_forces[i]);
}
for (i = 0; i < num_variables(); i++) {
previous_colvar_forces[i] = colvar_forces[i];
}
}
int colvarbias::end_of_step()
{
return COLVARS_OK;
}
int colvarbias::change_configuration(std::string const &conf)
{
cvm::error("Error: change_configuration() not implemented.\n",
COLVARS_NOT_IMPLEMENTED);
return COLVARS_NOT_IMPLEMENTED;
}
cvm::real colvarbias::energy_difference(std::string const &conf)
{
cvm::error("Error: energy_difference() not implemented.\n",
COLVARS_NOT_IMPLEMENTED);
return 0.0;
}
// So far, these are only implemented in colvarbias_abf
int colvarbias::bin_num()
{
cvm::error("Error: bin_num() not implemented.\n");
return COLVARS_NOT_IMPLEMENTED;
}
int colvarbias::current_bin()
{
cvm::error("Error: current_bin() not implemented.\n");
return COLVARS_NOT_IMPLEMENTED;
}
int colvarbias::bin_count(int bin_index)
{
cvm::error("Error: bin_count() not implemented.\n");
return COLVARS_NOT_IMPLEMENTED;
}
int colvarbias::replica_share()
{
cvm::error("Error: replica_share() not implemented.\n");
return COLVARS_NOT_IMPLEMENTED;
}
std::string const colvarbias::get_state_params() const
{
std::ostringstream os;
os << "step " << cvm::step_absolute() << "\n"
<< "name " << this->name << "\n";
return os.str();
}
int colvarbias::set_state_params(std::string const &conf)
{
std::string new_name = "";
if (colvarparse::get_keyval(conf, "name", new_name,
std::string(""), colvarparse::parse_silent) &&
(new_name != this->name)) {
cvm::error("Error: in the state file, the "
"\""+bias_type+"\" block has a different name, \""+new_name+
"\": different system?\n", INPUT_ERROR);
}
if (name.size() == 0) {
cvm::error("Error: \""+bias_type+"\" block within the restart file "
"has no identifiers.\n", INPUT_ERROR);
}
colvarparse::get_keyval(conf, "step", state_file_step,
cvm::step_absolute(), colvarparse::parse_silent);
return COLVARS_OK;
}
std::ostream & colvarbias::write_state(std::ostream &os)
{
if (cvm::debug()) {
cvm::log("Writing state file for bias \""+name+"\"\n");
}
os.setf(std::ios::scientific, std::ios::floatfield);
os.precision(cvm::cv_prec);
os << bias_type << " {\n"
<< " configuration {\n";
std::istringstream is(get_state_params());
std::string line;
while (std::getline(is, line)) {
os << " " << line << "\n";
}
os << " }\n";
write_state_data(os);
os << "}\n\n";
return os;
}
std::istream & colvarbias::read_state(std::istream &is)
{
size_t const start_pos = is.tellg();
std::string key, brace, conf;
if ( !(is >> key) || !(key == bias_type) ||
!(is >> brace) || !(brace == "{") ||
!(is >> colvarparse::read_block("configuration", conf)) ||
(set_state_params(conf) != COLVARS_OK) ) {
cvm::error("Error: in reading state configuration for \""+bias_type+"\" bias \""+
this->name+"\" at position "+
cvm::to_str(is.tellg())+" in stream.\n", INPUT_ERROR);
is.clear();
is.seekg(start_pos, std::ios::beg);
is.setstate(std::ios::failbit);
return is;
}
if (!read_state_data(is)) {
cvm::error("Error: in reading state data for \""+bias_type+"\" bias \""+
this->name+"\" at position "+
cvm::to_str(is.tellg())+" in stream.\n", INPUT_ERROR);
is.clear();
is.seekg(start_pos, std::ios::beg);
is.setstate(std::ios::failbit);
}
is >> brace;
if (brace != "}") {
cvm::error("Error: corrupt restart information for \""+bias_type+"\" bias \""+
this->name+"\": no matching brace at position "+
cvm::to_str(is.tellg())+" in stream.\n");
is.setstate(std::ios::failbit);
}
return is;
}
std::istream & colvarbias::read_state_data_key(std::istream &is, char const *key)
{
size_t const start_pos = is.tellg();
std::string key_in;
if ( !(is >> key_in) ||
!(key_in == to_lower_cppstr(std::string(key))) ) {
cvm::error("Error: in reading restart configuration for "+
bias_type+" bias \""+this->name+"\" at position "+
cvm::to_str(is.tellg())+" in stream.\n", INPUT_ERROR);
is.clear();
is.seekg(start_pos, std::ios::beg);
is.setstate(std::ios::failbit);
return is;
}
return is;
}
std::ostream & colvarbias::write_traj_label(std::ostream &os)
{
os << " ";
if (b_output_energy)
os << " E_"
<< cvm::wrap_string(this->name, cvm::en_width-2);
return os;
}
std::ostream & colvarbias::write_traj(std::ostream &os)
{
os << " ";
if (b_output_energy)
os << " "
<< std::setprecision(cvm::en_prec) << std::setw(cvm::en_width)
<< bias_energy;
return os;
}
colvarbias_ti::colvarbias_ti(char const *key)
: colvarbias(key)
{
provide(f_cvb_calc_ti_samples);
ti_avg_forces = NULL;
ti_count = NULL;
}
colvarbias_ti::~colvarbias_ti()
{
colvarbias_ti::clear_state_data();
}
int colvarbias_ti::clear_state_data()
{
if (ti_avg_forces != NULL) {
delete ti_avg_forces;
ti_avg_forces = NULL;
}
if (ti_count != NULL) {
delete ti_count;
ti_count = NULL;
}
return COLVARS_OK;
}
int colvarbias_ti::init(std::string const &conf)
{
int error_code = COLVARS_OK;
get_keyval_feature(this, conf, "writeTISamples",
f_cvb_write_ti_samples,
is_enabled(f_cvb_write_ti_samples));
get_keyval_feature(this, conf, "writeTIPMF",
f_cvb_write_ti_pmf,
is_enabled(f_cvb_write_ti_pmf));
if ((num_variables() > 1) && is_enabled(f_cvb_write_ti_pmf)) {
return cvm::error("Error: only 1-dimensional PMFs can be written "
"on the fly.\n"
"Consider using writeTISamples instead and "
"post-processing the sampled free-energy gradients.\n",
COLVARS_NOT_IMPLEMENTED);
} else {
error_code |= init_grids();
}
if (is_enabled(f_cvb_write_ti_pmf)) {
enable(f_cvb_write_ti_samples);
}
if (is_enabled(f_cvb_calc_ti_samples)) {
std::vector<std::string> const time_biases =
cvm::main()->time_dependent_biases();
if (time_biases.size() > 0) {
if ((time_biases.size() > 1) || (time_biases[0] != this->name)) {
for (size_t i = 0; i < num_variables(); i++) {
if (! variables(i)->is_enabled(f_cv_subtract_applied_force)) {
return cvm::error("Error: cannot collect TI samples while other "
"time-dependent biases are active and not all "
"variables have subtractAppliedForces on.\n",
INPUT_ERROR);
}
}
}
}
}
return error_code;
}
int colvarbias_ti::init_grids()
{
if (is_enabled(f_cvb_calc_ti_samples)) {
if (ti_avg_forces == NULL) {
ti_bin.resize(num_variables());
ti_system_forces.resize(num_variables());
for (size_t icv = 0; icv < num_variables(); icv++) {
ti_system_forces[icv].type(variables(icv)->value());
ti_system_forces[icv].is_derivative();
ti_system_forces[icv].reset();
}
ti_avg_forces = new colvar_grid_gradient(colvars);
ti_count = new colvar_grid_count(colvars);
ti_avg_forces->samples = ti_count;
ti_count->has_parent_data = true;
}
}
return COLVARS_OK;
}
int colvarbias_ti::update()
{
return update_system_forces(NULL);
}
int colvarbias_ti::update_system_forces(std::vector<colvarvalue> const
*subtract_forces)
{
if (! is_enabled(f_cvb_calc_ti_samples)) {
return COLVARS_OK;
}
has_data = true;
if (cvm::debug()) {
cvm::log("Updating system forces for bias "+this->name+"\n");
}
colvarproxy *proxy = cvm::main()->proxy;
size_t i;
if (proxy->total_forces_same_step()) {
for (i = 0; i < num_variables(); i++) {
ti_bin[i] = ti_avg_forces->current_bin_scalar(i);
}
}
// Collect total colvar forces
if ((cvm::step_relative() > 0) || proxy->total_forces_same_step()) {
if (ti_avg_forces->index_ok(ti_bin)) {
for (i = 0; i < num_variables(); i++) {
if (variables(i)->is_enabled(f_cv_subtract_applied_force)) {
// this colvar is already subtracting all applied forces
ti_system_forces[i] = variables(i)->total_force();
} else {
ti_system_forces[i] = variables(i)->total_force() -
((subtract_forces != NULL) ?
(*subtract_forces)[i] : previous_colvar_forces[i]);
}
}
ti_avg_forces->acc_value(ti_bin, ti_system_forces);
}
}
if (!proxy->total_forces_same_step()) {
// Set the index for use in the next iteration, when total forces come in
for (i = 0; i < num_variables(); i++) {
ti_bin[i] = ti_avg_forces->current_bin_scalar(i);
}
}
return COLVARS_OK;
}
std::string const colvarbias_ti::get_state_params() const
{
return std::string("");
}
int colvarbias_ti::set_state_params(std::string const &state_conf)
{
return COLVARS_OK;
}
std::ostream & colvarbias_ti::write_state_data(std::ostream &os)
{
if (! is_enabled(f_cvb_calc_ti_samples)) {
return os;
}
os << "\nhistogram\n";
ti_count->write_raw(os);
os << "\nsystem_forces\n";
ti_avg_forces->write_raw(os);
return os;
}
std::istream & colvarbias_ti::read_state_data(std::istream &is)
{
if (! is_enabled(f_cvb_calc_ti_samples)) {
return is;
}
if (cvm::debug()) {
cvm::log("Reading state data for the TI estimator.\n");
}
if (! read_state_data_key(is, "histogram")) {
return is;
}
if (! ti_count->read_raw(is)) {
return is;
}
if (! read_state_data_key(is, "system_forces")) {
return is;
}
if (! ti_avg_forces->read_raw(is)) {
return is;
}
if (cvm::debug()) {
cvm::log("Done reading state data for the TI estimator.\n");
}
return is;
}
int colvarbias_ti::write_output_files()
{
if (!has_data) {
// nothing to write
return COLVARS_OK;
}
std::string const ti_output_prefix = cvm::output_prefix()+"."+this->name;
std::ostream *os = NULL;
if (is_enabled(f_cvb_write_ti_samples)) {
std::string const ti_count_file_name(ti_output_prefix+".ti.count");
os = cvm::proxy->output_stream(ti_count_file_name);
if (os) {
ti_count->write_multicol(*os);
cvm::proxy->close_output_stream(ti_count_file_name);
}
std::string const ti_grad_file_name(ti_output_prefix+".ti.grad");
os = cvm::proxy->output_stream(ti_grad_file_name);
if (os) {
ti_avg_forces->write_multicol(*os);
cvm::proxy->close_output_stream(ti_grad_file_name);
}
}
if (is_enabled(f_cvb_write_ti_pmf)) {
std::string const pmf_file_name(ti_output_prefix+".ti.pmf");
cvm::log("Writing TI PMF to file \""+pmf_file_name+"\".\n");
os = cvm::proxy->output_stream(pmf_file_name);
if (os) {
// get the FE gradient
ti_avg_forces->multiply_constant(-1.0);
ti_avg_forces->write_1D_integral(*os);
ti_avg_forces->multiply_constant(-1.0);
cvm::proxy->close_output_stream(pmf_file_name);
}
}
return COLVARS_OK;
}
// Static members
std::vector<colvardeps::feature *> colvarbias::cvb_features;
|