#include <colvargrid.h>
Inheritance diagram for colvar_grid< T >:

Public Member Functions | |
| size_t | number_of_colvars () const |
| Return the number of colvars. | |
| size_t | number_of_points (int const icv=-1) const |
| std::vector< int > const & | sizes () const |
| Get the sizes in each direction. | |
| void | set_sizes (std::vector< int > const &new_sizes) |
| Set the sizes in each direction. | |
| size_t | multiplicity () const |
| Return the multiplicity of the type used. | |
| void | create (std::vector< int > const &nx_i, T const &t=T(), size_t const &mult_i=1) |
| Allocate data (allow initialization also after construction). | |
| void | create () |
| Allocate data (allow initialization also after construction). | |
| void | reset (T const &t=T()) |
| Reset data (in case the grid is being reused). | |
| colvar_grid () | |
| Default constructor. | |
| virtual | ~colvar_grid () |
| Destructor. | |
| colvar_grid (colvar_grid< T > const &g) | |
| "Almost copy-constructor": only copies configuration parameters from another grid, but doesn't reallocate stuff; create() must be called after that; | |
| colvar_grid (std::vector< int > const &nx_i, T const &t=T(), size_t const &mult_i=1) | |
| Constructor from explicit grid sizes. | |
| colvar_grid (std::vector< colvar * > const &colvars, T const &t=T(), size_t const &mult_i=1, bool margin=false) | |
| Constructor from a vector of colvars. | |
| void | wrap (std::vector< int > &ix) const |
| int | current_bin_scalar (int const i) const |
| Report the bin corresponding to the current value of variable i. | |
| int | value_to_bin_scalar (colvarvalue const &value, const int i) const |
| Use the lower boundary and the width to report which bin the provided value is in. | |
| int | value_to_bin_scalar (colvarvalue const &value, colvarvalue const &new_offset, cvm::real const &new_width) const |
| Same as the standard version, but uses another grid definition. | |
| colvarvalue | bin_to_value_scalar (int const &i_bin, int const i) const |
| Use the two boundaries and the width to report the central value corresponding to a bin index. | |
| colvarvalue | bin_to_value_scalar (int const &i_bin, colvarvalue const &new_offset, cvm::real const &new_width) const |
| Same as the standard version, but uses different parameters. | |
| void | set_value (std::vector< int > const &ix, T const &t, size_t const &imult=0) |
| Set the value at the point with index ix. | |
| T const & | value (std::vector< int > const &ix, size_t const &imult=0) const |
| Get the binned value indexed by ix, or the first of them if the multiplicity is larger than 1. | |
| void | add_constant (T const &t) |
| Add a constant to all elements (fast loop). | |
| void | multiply_constant (cvm::real const &a) |
| Multiply all elements by a scalar constant (fast loop). | |
| std::vector< int > const | get_colvars_index (std::vector< colvarvalue > const &values) const |
| Get the bin indices corresponding to the provided values of the colvars. | |
| std::vector< int > const | get_colvars_index () const |
| Get the bin indices corresponding to the current values of the colvars. | |
| cvm::real | bin_distance_from_boundaries (std::vector< colvarvalue > const &values) |
| Get the minimal distance (in number of bins) from the boundaries; a negative number is returned if the given point is off-grid. | |
| void | map_grid (colvar_grid< T > const &other_grid) |
| Add data from another grid of the same type. | |
| void | add_grid (colvar_grid< T > const &other_grid, cvm::real scale_factor=1.0) |
| Add data from another grid of the same type, AND identical definition (boundaries, widths). | |
| virtual T | value_output (std::vector< int > const &ix, size_t const &imult=0) |
| Return the value suitable for output purposes (so that it may be rescaled or manipulated without changing it permanently). | |
| virtual void | value_input (std::vector< int > const &ix, T const &t, size_t const &imult=0, bool add=false) |
| Get the value from a formatted output and transform it into the internal representation (the two may be different, e.g. when using colvar_grid_count). | |
| std::vector< int > const | new_index () const |
| Get the index corresponding to the "first" bin, to be used as the initial value for an index in looping. | |
| bool | index_ok (std::vector< int > const &ix) const |
| Check that the index is within range in each of the dimensions. | |
| void | incr (std::vector< int > &ix) const |
| Increment the index, in a way that will make it loop over the whole nd-dimensional array. | |
| std::ostream & | write_params (std::ostream &os) |
| Write the grid parameters (number of colvars, boundaries, width and number of points). | |
| bool | parse_params (std::string const &conf) |
| void | check_consistency () |
| Check that the grid information inside (boundaries, widths, ...) is consistent with the current setting of the colvars. | |
| void | check_consistency (colvar_grid< T > const &other_grid) |
| Check that the grid information inside (boundaries, widths, ...) is consistent with the one of another grid. | |
| std::ostream & | write_raw (std::ostream &os, size_t const buf_size=3) |
| Write the grid data without labels, as they are represented in memory. | |
| std::istream & | read_raw (std::istream &is) |
| Read data written by colvar_grid::write_raw(). | |
| void | read_raw_error () |
| To be called after colvar_grid::read_raw() returns an error. | |
| void | write_multicol (std::ostream &os) |
| Write the grid in a format which is both human readable and suitable for visualization e.g. with gnuplot. | |
| std::istream & | read_multicol (std::istream &is, bool add=false) |
| Read a grid written by colvar_grid::write_multicol() Adding data if add is true, replacing if false. | |
Public Attributes | |
| std::vector< colvarvalue > | lower_boundaries |
| Lower boundaries of the colvars in this grid. | |
| std::vector< colvarvalue > | upper_boundaries |
| Upper boundaries of the colvars in this grid. | |
| std::vector< bool > | periodic |
| Whether some colvars are periodic. | |
| std::vector< cvm::real > | widths |
| Widths of the colvars in this grid. | |
| bool | has_parent_data |
| True if this is a count grid related to another grid of data. | |
| bool | has_data |
| Whether this grid has been filled with data or is still empty. | |
Protected Member Functions | |
| size_t | address (std::vector< int > const &ix) const |
| Get the low-level index corresponding to an index. | |
Protected Attributes | |
| size_t | nd |
| Number of dimensions. | |
| std::vector< int > | nx |
| Number of points along each dimension. | |
| std::vector< int > | nxc |
| Cumulative number of points along each dimension. | |
| size_t | mult |
| Multiplicity of each datum (allow the binning of non-scalar types). | |
| size_t | nt |
| Total number of grid points. | |
| std::vector< T > | data |
| Low-level array of values. | |
| std::vector< size_t > | new_data |
| Newly read data (used for count grids, when adding several grids read from disk). | |
| std::vector< colvar * > | cv |
| Colvars collected in this grid. | |
| T | The data type |
Definition at line 19 of file colvargrid.h.
|
|||||||||
|
Default constructor.
Definition at line 154 of file colvargrid.h. 00154 : has_data (false) 00155 { 00156 save_delimiters = false; 00157 nd = nt = 0; 00158 }
|
|
|||||||||
|
Destructor.
Definition at line 161 of file colvargrid.h. 00162 {}
|
|
||||||||||
|
"Almost copy-constructor": only copies configuration parameters from another grid, but doesn't reallocate stuff; create() must be called after that;
Definition at line 167 of file colvargrid.h. 00167 : has_data (false), 00168 nd (g.nd), 00169 mult (g.mult), 00170 cv (g.cv), 00171 lower_boundaries (g.lower_boundaries), 00172 upper_boundaries (g.upper_boundaries), 00173 periodic (g.periodic), 00174 widths (g.widths), 00175 data() 00176 { 00177 save_delimiters = false; 00178 }
|
|
||||||||||||||||||||
|
Constructor from explicit grid sizes.
Definition at line 184 of file colvargrid.h. 00186 : has_data (false) 00187 { 00188 save_delimiters = false; 00189 this->create (nx_i, t, mult_i); 00190 }
|
|
||||||||||||||||||||||||
|
Constructor from a vector of colvars.
Definition at line 193 of file colvargrid.h. 00197 : cv (colvars), has_data (false) 00198 { 00199 save_delimiters = false; 00200 00201 std::vector<int> nx_i; 00202 00203 if (cvm::debug()) 00204 cvm::log ("Allocating a grid for "+cvm::to_str (colvars.size())+ 00205 " collective variables.\n"); 00206 00207 for (size_t i = 0; i < cv.size(); i++) { 00208 00209 if (cv[i]->type() != colvarvalue::type_scalar) { 00210 cvm::fatal_error ("Colvar grids can only be automatically " 00211 "constructed for scalar variables. " 00212 "ABF and histogram can not be used; metadynamics " 00213 "can be used with useGrids disabled.\n"); 00214 } 00215 00216 if (cv[i]->width <= 0.0) { 00217 cvm::fatal_error ("Tried to initialize a grid on a " 00218 "variable with negative or zero width.\n"); 00219 } 00220 00221 if (!cv[i]->tasks[colvar::task_lower_boundary] || !cv[i]->tasks[colvar::task_upper_boundary]) { 00222 cvm::fatal_error ("Tried to initialize a grid on a " 00223 "variable with undefined boundaries.\n"); 00224 } 00225 00226 widths.push_back (cv[i]->width); 00227 periodic.push_back (cv[i]->periodic_boundaries()); 00228 00229 if (margin) { 00230 if (periodic[i]) { 00231 // Shift the grid by half the bin width (values at edges instead of center of bins) 00232 lower_boundaries.push_back (cv[i]->lower_boundary.real_value - 0.5 * widths[i]); 00233 upper_boundaries.push_back (cv[i]->upper_boundary.real_value - 0.5 * widths[i]); 00234 } else { 00235 // Make this grid larger by one bin width 00236 lower_boundaries.push_back (cv[i]->lower_boundary.real_value - 0.5 * widths[i]); 00237 upper_boundaries.push_back (cv[i]->upper_boundary.real_value + 0.5 * widths[i]); 00238 } 00239 } else { 00240 lower_boundaries.push_back (cv[i]->lower_boundary); 00241 upper_boundaries.push_back (cv[i]->upper_boundary); 00242 } 00243 00244 00245 { 00246 cvm::real nbins = ( upper_boundaries[i].real_value - 00247 lower_boundaries[i].real_value ) / widths[i]; 00248 int nbins_round = (int)(nbins+0.5); 00249 00250 if (std::fabs (nbins - cvm::real (nbins_round)) > 1.0E-10) { 00251 cvm::log ("Warning: grid interval ("+ 00252 cvm::to_str (lower_boundaries[i], cvm::cv_width, cvm::cv_prec)+" - "+ 00253 cvm::to_str (upper_boundaries[i], cvm::cv_width, cvm::cv_prec)+ 00254 ") is not commensurate to its bin width ("+ 00255 cvm::to_str (widths[i], cvm::cv_width, cvm::cv_prec)+").\n"); 00256 upper_boundaries[i].real_value = lower_boundaries[i].real_value + 00257 (nbins_round * widths[i]); 00258 } 00259 00260 if (cvm::debug()) 00261 cvm::log ("Number of points is "+cvm::to_str ((int) nbins_round)+ 00262 " for the colvar no. "+cvm::to_str (i+1)+".\n"); 00263 00264 nx_i.push_back (nbins_round); 00265 } 00266 00267 } 00268 00269 create (nx_i, t, mult_i); 00270 }
|
|
||||||||||
|
Add a constant to all elements (fast loop).
Definition at line 344 of file colvargrid.h. Referenced by colvarbias_meta::write_pmf().
|
|
||||||||||||||||
|
Add data from another grid of the same type, AND identical definition (boundaries, widths).
Definition at line 454 of file colvargrid.h. Referenced by colvarbias_meta::write_pmf(). 00456 {
00457 if (other_grid.multiplicity() != this->multiplicity())
00458 cvm::fatal_error ("Error: trying to sum togetehr two grids with values of "
00459 "different multiplicity.\n");
00460 if (scale_factor != 1.0)
00461 for (size_t i = 0; i < data.size(); i++) {
00462 data[i] += scale_factor * other_grid.data[i];
00463 }
00464 else
00465 // skip multiplication if possible
00466 for (size_t i = 0; i < data.size(); i++) {
00467 data[i] += other_grid.data[i];
00468 }
00469 has_data = true;
00470 }
|
|
||||||||||
|
Get the low-level index corresponding to an index.
Definition at line 49 of file colvargrid.h. Referenced by colvar_grid< size_t >::set_value(), colvar_grid< size_t >::value(), and colvar_grid< size_t >::value_input(). 00050 {
00051 size_t addr = 0;
00052 for (size_t i = 0; i < nd; i++) {
00053 addr += ix[i]*nxc[i];
00054 if (cvm::debug()) {
00055 if (ix[i] >= nx[i])
00056 cvm::fatal_error ("Error: exceeding bounds in colvar_grid.\n");
00057 }
00058 }
00059 return addr;
00060 }
|
|
||||||||||
|
Get the minimal distance (in number of bins) from the boundaries; a negative number is returned if the given point is off-grid.
Definition at line 384 of file colvargrid.h. Referenced by colvarbias_meta::create_hill(), colvarbias_meta::read_hill(), and colvarbias_meta::recount_hills_off_grid(). 00385 {
00386 cvm::real minimum = 1.0E+16;
00387 for (size_t i = 0; i < nd; i++) {
00388
00389 if (periodic[i]) continue;
00390
00391 cvm::real dl = std::sqrt (cv[i]->dist2 (values[i], lower_boundaries[i])) / widths[i];
00392 cvm::real du = std::sqrt (cv[i]->dist2 (values[i], upper_boundaries[i])) / widths[i];
00393
00394 if (values[i].real_value < lower_boundaries[i])
00395 dl *= -1.0;
00396 if (values[i].real_value > upper_boundaries[i])
00397 du *= -1.0;
00398
00399 if (dl < minimum)
00400 minimum = dl;
00401 if (du < minimum)
00402 minimum = du;
00403 }
00404
00405 return minimum;
00406 }
|
|
||||||||||||||||||||
|
Same as the standard version, but uses different parameters.
Definition at line 317 of file colvargrid.h. 00320 {
00321 return new_offset.real_value + new_width * (0.5 + i_bin);
00322 }
|
|
||||||||||||||||
|
Use the two boundaries and the width to report the central value corresponding to a bin index.
Definition at line 311 of file colvargrid.h. Referenced by colvarbias_meta::project_hills(), and colvar_grid< size_t >::write_multicol(). 00312 {
00313 return lower_boundaries[i].real_value + widths[i] * (0.5 + i_bin);
00314 }
|
|
||||||||||
|
Check that the grid information inside (boundaries, widths, ...) is consistent with the one of another grid.
Definition at line 640 of file colvargrid.h. 00641 {
00642 for (size_t i = 0; i < nd; i++) {
00643 // we skip dist2(), because periodicities and the like should
00644 // matter: boundaries should be EXACTLY the same (otherwise,
00645 // map_grid() should be used)
00646 if ( (std::fabs (other_grid.lower_boundaries[i] -
00647 lower_boundaries[i]) > 1.0E-10) ||
00648 (std::fabs (other_grid.upper_boundaries[i] -
00649 upper_boundaries[i]) > 1.0E-10) ||
00650 (std::fabs (other_grid.widths[i] -
00651 widths[i]) > 1.0E-10) ||
00652 (data.size() != other_grid.data.size()) ) {
00653 cvm::fatal_error ("Error: inconsistency between "
00654 "two grids that are supposed to be equal, "
00655 "aside from the data stored.\n");
00656 }
00657 }
00658 }
|
|
|||||||||
|
Check that the grid information inside (boundaries, widths, ...) is consistent with the current setting of the colvars.
Definition at line 622 of file colvargrid.h. 00623 {
00624 for (size_t i = 0; i < nd; i++) {
00625 if ( (std::sqrt (cv[i]->dist2 (cv[i]->lower_boundary,
00626 lower_boundaries[i])) > 1.0E-10) ||
00627 (std::sqrt (cv[i]->dist2 (cv[i]->upper_boundary,
00628 upper_boundaries[i])) > 1.0E-10) ||
00629 (std::sqrt (cv[i]->dist2 (cv[i]->width,
00630 widths[i])) > 1.0E-10) ) {
00631 cvm::fatal_error ("Error: restart information for a grid is "
00632 "inconsistent with that of its colvars.\n");
00633 }
00634 }
00635 }
|
|
|||||||||
|
Allocate data (allow initialization also after construction).
Definition at line 141 of file colvargrid.h. Referenced by colvar_grid< size_t >::colvar_grid(), colvar_grid< size_t >::create(), and colvar_grid< size_t >::parse_params(). 00142 {
00143 create (this->nx, T(), this->mult);
00144 }
|
|
||||||||||||||||||||
|
Allocate data (allow initialization also after construction).
Definition at line 118 of file colvargrid.h. Referenced by colvarbias_meta::update(). 00121 {
00122 mult = mult_i;
00123 nd = nx_i.size();
00124 nxc.resize (nd);
00125 nx = nx_i;
00126
00127 nt = mult;
00128 for (int i = nd-1; i >= 0; i--) {
00129 if (nx_i[i] <= 0)
00130 cvm::fatal_error ("Error: providing an invalid number of points, "+
00131 cvm::to_str (nx_i[i])+".\n");
00132 nxc[i] = nt;
00133 nt *= nx[i];
00134 }
00135
00136 data.reserve (nt);
00137 data.assign (nt, t);
00138 }
|
|
||||||||||
|
Report the bin corresponding to the current value of variable i.
Definition at line 289 of file colvargrid.h. Referenced by colvar_grid< size_t >::get_colvars_index(), colvarbias_histogram::update(), and colvarbias_abf::update(). 00290 {
00291 return value_to_bin_scalar (cv[i]->value(), i);
00292 }
|
|
|||||||||
|
Get the bin indices corresponding to the current values of the colvars.
Definition at line 372 of file colvargrid.h. 00373 {
00374 std::vector<int> index = new_index();
00375 for (size_t i = 0; i < nd; i++) {
00376 index[i] = current_bin_scalar (i);
00377 }
00378 return index;
00379 }
|
|
||||||||||
|
Get the bin indices corresponding to the provided values of the colvars.
Definition at line 361 of file colvargrid.h. Referenced by colvarbias_meta::update(). 00362 {
00363 std::vector<int> index = new_index();
00364 for (size_t i = 0; i < nd; i++) {
00365 index[i] = value_to_bin_scalar (values[i], i);
00366 }
00367 return index;
00368 }
|
|
||||||||||
|
Increment the index, in a way that will make it loop over the whole nd-dimensional array.
Definition at line 521 of file colvargrid.h. Referenced by colvar_grid< size_t >::map_grid(), colvarbias_meta::project_hills(), colvar_grid< size_t >::read_multicol(), colvar_grid< size_t >::read_raw(), colvar_grid< size_t >::write_multicol(), and colvar_grid< size_t >::write_raw(). 00522 {
00523 for (int i = ix.size()-1; i >= 0; i--) {
00524
00525 ix[i]++;
00526
00527 if (ix[i] >= nx[i]) {
00528
00529 if (i > 0) {
00530 ix[i] = 0;
00531 continue;
00532 } else {
00533 // this is the last iteration, a non-valid index is being
00534 // set for the outer index, which will be caught by
00535 // index_ok()
00536 ix[0] = nx[0];
00537 return;
00538 }
00539 } else {
00540 return;
00541 }
00542 }
00543 }
|
|
||||||||||
|
Check that the index is within range in each of the dimensions.
Definition at line 510 of file colvargrid.h. Referenced by colvar_grid< size_t >::map_grid(), colvarbias_meta::project_hills(), colvar_grid< size_t >::read_multicol(), colvar_grid< size_t >::read_raw(), colvarbias_meta::update(), colvarbias_histogram::update(), colvarbias_abf::update(), colvar_grid< size_t >::write_multicol(), and colvar_grid< size_t >::write_raw(). 00511 {
00512 for (size_t i = 0; i < nd; i++) {
00513 if ( (ix[i] < 0) || (ix[i] >= int (nx[i])) )
00514 return false;
00515 }
00516 return true;
00517 }
|
|
||||||||||
|
Add data from another grid of the same type. Note: this function maps other_grid inside this one regardless of whether it fits or not. Definition at line 413 of file colvargrid.h. Referenced by colvarbias_meta::read_restart(), and colvarbias_meta::update(). 00414 {
00415 if (other_grid.multiplicity() != this->multiplicity())
00416 cvm::fatal_error ("Error: trying to merge two grids with values of "
00417 "different multiplicity.\n");
00418
00419 std::vector<colvarvalue> const &gb = this->lower_boundaries;
00420 std::vector<cvm::real> const &gw = this->widths;
00421 std::vector<colvarvalue> const &ogb = other_grid.lower_boundaries;
00422 std::vector<cvm::real> const &ogw = other_grid.widths;
00423
00424 std::vector<int> ix = this->new_index();
00425 std::vector<int> oix = other_grid.new_index();
00426
00427 if (cvm::debug())
00428 cvm::log ("Remapping grid...\n");
00429 for ( ; this->index_ok (ix); this->incr (ix)) {
00430
00431 for (size_t i = 0; i < nd; i++) {
00432 oix[i] =
00433 value_to_bin_scalar (bin_to_value_scalar (ix[i], gb[i], gw[i]),
00434 ogb[i],
00435 ogw[i]);
00436 }
00437
00438 if (! other_grid.index_ok (oix)) {
00439 continue;
00440 }
00441
00442 for (size_t im = 0; im < mult; im++) {
00443 this->set_value (ix, other_grid.value (oix, im), im);
00444 }
00445 }
00446
00447 has_data = true;
00448 if (cvm::debug())
00449 cvm::log ("Remapping done.\n");
00450 }
|
|
|||||||||
|
Return the multiplicity of the type used.
Definition at line 112 of file colvargrid.h. Referenced by colvar_grid< size_t >::add_grid(), and colvar_grid< size_t >::map_grid(). 00113 {
00114 return mult;
00115 }
|
|
||||||||||
|
Multiply all elements by a scalar constant (fast loop).
Definition at line 352 of file colvargrid.h. Referenced by colvarbias_meta::write_pmf(). 00353 {
00354 for (size_t i = 0; i < nt; i++)
00355 data[i] *= a;
00356 }
|
|
|||||||||
|
Get the index corresponding to the "first" bin, to be used as the initial value for an index in looping.
Definition at line 503 of file colvargrid.h. Referenced by colvar_grid< size_t >::get_colvars_index(), colvar_grid< size_t >::map_grid(), colvarbias_meta::project_hills(), colvar_grid< size_t >::read_multicol(), colvar_grid< size_t >::read_raw(), colvar_grid< size_t >::write_multicol(), and colvar_grid< size_t >::write_raw(). 00504 {
00505 return std::vector<int> (nd, 0);
00506 }
|
|
|||||||||
|
Return the number of colvars.
Definition at line 83 of file colvargrid.h. 00084 {
00085 return nd;
00086 }
|
|
||||||||||
|
Return the number of points in the i-th direction, if provided, or the total number Definition at line 90 of file colvargrid.h. 00091 {
00092 if (icv < 0) {
00093 return nt;
00094 } else {
00095 return nx[icv];
00096 }
00097 }
|
|
||||||||||
|
Definition at line 575 of file colvargrid.h. 00576 {
00577 std::vector<int> old_nx = nx;
00578 std::vector<colvarvalue> old_lb = lower_boundaries;
00579
00580 {
00581 size_t nd_in = 0;
00582 colvarparse::get_keyval (conf, "n_colvars", nd_in, nd, colvarparse::parse_silent);
00583 if (nd_in != nd)
00584 cvm::fatal_error ("Error: trying to read data for a grid "
00585 "that contains a different number of colvars ("+
00586 cvm::to_str (nd_in)+") than the grid defined "
00587 "in the configuration file ("+cvm::to_str (nd)+
00588 ").\n");
00589 }
00590
00591 colvarparse::get_keyval (conf, "lower_boundaries",
00592 lower_boundaries, lower_boundaries, colvarparse::parse_silent);
00593
00594 colvarparse::get_keyval (conf, "upper_boundaries",
00595 upper_boundaries, upper_boundaries, colvarparse::parse_silent);
00596
00597 colvarparse::get_keyval (conf, "widths", widths, widths, colvarparse::parse_silent);
00598
00599 colvarparse::get_keyval (conf, "sizes", nx, nx, colvarparse::parse_silent);
00600
00601 bool new_params = false;
00602 for (size_t i = 0; i < nd; i++) {
00603 if ( (old_nx[i] != nx[i]) ||
00604 (std::sqrt (cv[i]->dist2 (old_lb[i],
00605 lower_boundaries[i])) > 1.0E-10) ) {
00606 new_params = true;
00607 }
00608 }
00609
00610 // reallocate the array in case the grid params have just changed
00611 if (new_params) {
00612 data.resize (0);
00613 this->create (nx, T(), mult);
00614 }
00615
00616 return true;
00617 }
|
|
||||||||||||||||
|
Read a grid written by colvar_grid::write_multicol() Adding data if add is true, replacing if false.
Definition at line 760 of file colvargrid.h. 00761 {
00762 // Data in the header: nColvars, then for each
00763 // xiMin, dXi, nPoints, periodic
00764
00765 std::string hash;
00766 cvm::real lower, width, x;
00767 size_t n, periodic;
00768 bool remap;
00769 std::vector<T> new_value;
00770 std::vector<int> nx_read;
00771 std::vector<int> bin;
00772
00773 if ( cv.size() != nd ) {
00774 cvm::fatal_error ("Cannot read grid file: missing reference to colvars.");
00775 }
00776
00777 if ( !(is >> hash) || (hash != "#") ) {
00778 cvm::fatal_error ("Error reading grid at position "+
00779 cvm::to_str (is.tellg())+" in stream (read \"" + hash + "\")\n");
00780 }
00781
00782 is >> n;
00783 if ( n != nd ) {
00784 cvm::fatal_error ("Error reading grid: wrong number of collective variables.\n");
00785 }
00786
00787 nx_read.resize (n);
00788 bin.resize (n);
00789 new_value.resize (mult);
00790
00791 if (this->has_parent_data && add) {
00792 new_data.resize (data.size());
00793 }
00794
00795 remap = false;
00796 for (size_t i = 0; i < nd; i++ ) {
00797 if ( !(is >> hash) || (hash != "#") ) {
00798 cvm::fatal_error ("Error reading grid at position "+
00799 cvm::to_str (is.tellg())+" in stream (read \"" + hash + "\")\n");
00800 }
00801
00802 is >> lower >> width >> nx_read[i] >> periodic;
00803
00804
00805 if ( (std::fabs (lower - lower_boundaries[i].real_value) > 1.0e-10) ||
00806 (std::fabs (width - widths[i] ) > 1.0e-10) ||
00807 (nx_read[i] != nx[i]) ) {
00808 cvm::log ("Warning: reading from different grid definition (colvar "
00809 + cvm::to_str (i+1) + "); remapping data on new grid.\n");
00810 remap = true;
00811 }
00812 }
00813
00814 if ( remap ) {
00815 // re-grid data
00816 while (is.good()) {
00817 bool end_of_file = false;
00818
00819 for (size_t i = 0; i < nd; i++ ) {
00820 if ( !(is >> x) ) end_of_file = true;
00821 bin[i] = value_to_bin_scalar (x, i);
00822 }
00823 if (end_of_file) break;
00824
00825 for (size_t imult = 0; imult < mult; imult++) {
00826 is >> new_value[imult];
00827 }
00828
00829 if ( index_ok(bin) ) {
00830 for (size_t imult = 0; imult < mult; imult++) {
00831 value_input (bin, new_value[imult], imult, add);
00832 }
00833 }
00834 }
00835 } else {
00836 // do not re-grid the data but assume the same grid is used
00837 for (std::vector<int> ix = new_index(); index_ok (ix); incr (ix) ) {
00838 for (size_t i = 0; i < nd; i++ ) {
00839 is >> x;
00840 }
00841 for (size_t imult = 0; imult < mult; imult++) {
00842 is >> new_value[imult];
00843 value_input (ix, new_value[imult], imult, add);
00844 }
00845 }
00846 }
00847 has_data = true;
00848 return is;
00849 }
|
|
||||||||||
|
Read data written by colvar_grid::write_raw().
Definition at line 689 of file colvargrid.h. 00690 {
00691 size_t const start_pos = is.tellg();
00692
00693 for (std::vector<int> ix = new_index(); index_ok (ix); incr (ix)) {
00694 for (size_t imult = 0; imult < mult; imult++) {
00695 T new_value;
00696 if (is >> new_value) {
00697 value_input (ix, new_value, imult);
00698 } else {
00699 is.clear();
00700 is.seekg (start_pos, std::ios::beg);
00701 is.setstate (std::ios::failbit);
00702 return is;
00703 }
00704 }
00705 }
00706
00707 has_data = true;
00708 return is;
00709 }
|
|
|||||||||
|
To be called after colvar_grid::read_raw() returns an error.
Definition at line 712 of file colvargrid.h. 00713 {
00714 cvm::fatal_error ("Error: failed to read all of the grid points from file. Possible explanations: grid parameters in the configuration (lowerBoundary, upperBoundary, width) are different from those in the file, or the file is corrupt/incomplete.\n");
00715 }
|
|
||||||||||
|
Reset data (in case the grid is being reused).
Definition at line 147 of file colvargrid.h. Referenced by colvarbias_meta::project_hills(), colvarbias_meta::update(), colvarbias_abf::update(), and colvarbias_meta::write_pmf(). 00148 {
00149 data.assign (nt, t);
00150 }
|
|
||||||||||
|
Set the sizes in each direction.
Definition at line 106 of file colvargrid.h. 00107 {
00108 nx = new_sizes;
00109 }
|
|
||||||||||||||||||||
|
Set the value at the point with index ix.
Definition at line 325 of file colvargrid.h. Referenced by colvar_grid< size_t >::map_grid().
|
|
|||||||||
|
Get the sizes in each direction.
Definition at line 100 of file colvargrid.h. Referenced by colvarbias_meta::update(). 00101 {
00102 return nx;
00103 }
|
|
||||||||||||||||
|
Get the binned value indexed by ix, or the first of them if the multiplicity is larger than 1.
Definition at line 336 of file colvargrid.h. Referenced by colvar_grid_gradient::average(), colvarbias_meta::calc_hills(), colvar_grid< size_t >::current_bin_scalar(), colvar_grid< size_t >::map_grid(), colvarbias_meta::update(), colvarbias_abf::update(), colvar_grid_gradient::value_input(), colvar_grid_scalar::value_input(), colvar_grid_gradient::value_output(), colvar_grid_scalar::value_output(), colvar_grid< size_t >::value_output(), and colvar_grid_gradient::write_1D_integral().
|
|
||||||||||||||||||||||||
|
Get the value from a formatted output and transform it into the internal representation (the two may be different, e.g. when using colvar_grid_count).
Reimplemented in colvar_grid_count, colvar_grid_scalar, and colvar_grid_gradient. Definition at line 483 of file colvargrid.h. Referenced by colvar_grid< size_t >::read_multicol(), and colvar_grid< size_t >::read_raw(). 00487 {
00488 if ( add )
00489 data[address (ix) + imult] += t;
00490 else
00491 data[address (ix) + imult] = t;
00492 has_data = true;
00493 }
|
|
||||||||||||||||
|
Return the value suitable for output purposes (so that it may be rescaled or manipulated without changing it permanently).
Reimplemented in colvar_grid_scalar, and colvar_grid_gradient. Definition at line 474 of file colvargrid.h. Referenced by colvar_grid< size_t >::write_multicol(), and colvar_grid< size_t >::write_raw(). 00476 {
00477 return value (ix, imult);
00478 }
|
|
||||||||||||||||||||
|
Same as the standard version, but uses another grid definition.
Definition at line 302 of file colvargrid.h. 00305 {
00306 return (int) std::floor ( (value.real_value - new_offset.real_value) / new_width );
00307 }
|
|
||||||||||||||||
|
Use the lower boundary and the width to report which bin the provided value is in.
Definition at line 296 of file colvargrid.h. Referenced by colvar_grid< size_t >::current_bin_scalar(), colvar_grid< size_t >::get_colvars_index(), colvar_grid< size_t >::map_grid(), and colvar_grid< size_t >::read_multicol(). 00297 {
00298 return (int) std::floor ( (value.real_value - lower_boundaries[i].real_value) / widths[i] );
00299 }
|
|
||||||||||
|
Wrap an index vector around periodic boundary conditions also checks validity of non-periodic indices Definition at line 275 of file colvargrid.h. 00276 {
00277 for (size_t i = 0; i < nd; i++) {
00278 if (periodic[i]) {
00279 ix[i] = (ix[i] + nx[i]) % nx[i]; //to ensure non-negative result
00280 } else {
00281 if (ix[i] < 0 || ix[i] >= nx[i])
00282 cvm::fatal_error ("Trying to wrap illegal index vector (non-PBC): "
00283 + cvm::to_str (ix));
00284 }
00285 }
00286 }
|
|
||||||||||
|
Write the grid in a format which is both human readable and suitable for visualization e.g. with gnuplot.
Definition at line 719 of file colvargrid.h. Referenced by colvarbias_histogram::update(), and colvarbias_meta::write_pmf(). 00720 {
00721 std::streamsize const w = os.width();
00722 std::streamsize const p = os.precision();
00723
00724 // Data in the header: nColvars, then for each
00725 // xiMin, dXi, nPoints, periodic
00726
00727 os << std::setw (2) << "# " << nd << "\n";
00728 for (size_t i = 0; i < nd; i++) {
00729 os << "# "
00730 << std::setw (10) << lower_boundaries[i]
00731 << std::setw (10) << widths[i]
00732 << std::setw (10) << nx[i] << " "
00733 << periodic[i] << "\n";
00734 }
00735
00736 for (std::vector<int> ix = new_index(); index_ok (ix); incr (ix) ) {
00737
00738 if (ix.back() == 0) {
00739 // if the last index is 0, add a new line to mark the new record
00740 os << "\n";
00741 }
00742
00743 for (size_t i = 0; i < nd; i++) {
00744 os << " "
00745 << std::setw (w) << std::setprecision (p)
00746 << bin_to_value_scalar (ix[i], i);
00747 }
00748 os << " ";
00749 for (size_t imult = 0; imult < mult; imult++) {
00750 os << " "
00751 << std::setw (w) << std::setprecision (p)
00752 << value_output (ix, imult);
00753 }
00754 os << "\n";
00755 }
00756 }
|
|
||||||||||
|
Write the grid parameters (number of colvars, boundaries, width and number of points).
Definition at line 546 of file colvargrid.h. 00547 {
00548 os << "grid_parameters {\n n_colvars " << nd << "\n";
00549
00550 os << " lower_boundaries ";
00551 for (size_t i = 0; i < nd; i++)
00552 os << " " << lower_boundaries[i];
00553 os << "\n";
00554
00555 os << " upper_boundaries ";
00556 for (size_t i = 0; i < nd; i++)
00557 os << " " << upper_boundaries[i];
00558 os << "\n";
00559
00560 os << " widths ";
00561 for (size_t i = 0; i < nd; i++)
00562 os << " " << widths[i];
00563 os << "\n";
00564
00565 os << " sizes ";
00566 for (size_t i = 0; i < nd; i++)
00567 os << " " << nx[i];
00568 os << "\n";
00569
00570 os << "}\n";
00571 return os;
00572 }
|
|
||||||||||||||||
|
Write the grid data without labels, as they are represented in memory.
Definition at line 664 of file colvargrid.h. 00666 {
00667 std::streamsize const w = os.width();
00668 std::streamsize const p = os.precision();
00669
00670 std::vector<int> ix = new_index();
00671 size_t count = 0;
00672 for ( ; index_ok (ix); incr (ix)) {
00673 for (size_t imult = 0; imult < mult; imult++) {
00674 os << " "
00675 << std::setw (w) << std::setprecision (p)
00676 << value_output (ix, imult);
00677 if (((++count) % buf_size) == 0)
00678 os << "\n";
00679 }
00680 }
00681 // write a final newline only if buffer is not empty
00682 if ((count % buf_size) != 0)
00683 os << "\n";
00684
00685 return os;
00686 }
|
|
|||||
|
Colvars collected in this grid.
Definition at line 46 of file colvargrid.h. |
|
|||||
|
Low-level array of values.
Definition at line 40 of file colvargrid.h. Referenced by colvar_grid< size_t >::add_grid(), and colvar_grid< size_t >::check_consistency(). |
|
|||||
|
Whether this grid has been filled with data or is still empty.
Definition at line 80 of file colvargrid.h. |
|
|||||
|
True if this is a count grid related to another grid of data.
Definition at line 77 of file colvargrid.h. Referenced by colvarbias_abf::colvarbias_abf(). |
|
|||||
|
Lower boundaries of the colvars in this grid.
Definition at line 65 of file colvargrid.h. Referenced by colvar_grid< size_t >::check_consistency(), colvar_grid< size_t >::map_grid(), and colvarbias_meta::update(). |
|
|||||
|
Multiplicity of each datum (allow the binning of non-scalar types).
Definition at line 34 of file colvargrid.h. |
|
|||||
|
Number of dimensions.
Definition at line 24 of file colvargrid.h. |
|
|||||
|
Newly read data (used for count grids, when adding several grids read from disk).
Definition at line 43 of file colvargrid.h. |
|
|||||
|
Total number of grid points.
Definition at line 37 of file colvargrid.h. |
|
|||||
|
Number of points along each dimension.
Definition at line 27 of file colvargrid.h. |
|
|||||
|
Cumulative number of points along each dimension.
Definition at line 30 of file colvargrid.h. |
|
|||||
|
Whether some colvars are periodic.
Definition at line 71 of file colvargrid.h. |
|
|||||
|
Upper boundaries of the colvars in this grid.
Definition at line 68 of file colvargrid.h. Referenced by colvar_grid< size_t >::check_consistency(), and colvarbias_meta::update(). |
|
|||||
|
Widths of the colvars in this grid.
Definition at line 74 of file colvargrid.h. Referenced by colvarbias_meta::calc_hills(), colvar_grid< size_t >::check_consistency(), and colvar_grid< size_t >::map_grid(). |
1.3.9.1