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

colvar_grid< T > Class Template Reference

Grid of values of a function of several collective variables. More...

#include <colvargrid.h>

Inheritance diagram for colvar_grid< T >:

colvarparse List of all members.

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).
 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 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.
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.
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 minimum distance (in number of bins) from the boundaries; give a negative number if the point is off-grid.
void map_grid (colvar_grid< T > const &other_grid)
 Add data from another grid of the same type.
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 (it may have been rescaled or manipulated).
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.
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 write_multicol (std::ostream &os)
 Write the grid in a format which is both human readable and suitable for visualization e.g. with gnuplot.
void 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< colvarvaluelower_boundaries
 Lower boundaries of the colvars in this grid.
std::vector< colvarvalueupper_boundaries
 Upper boundaries of the colvars in this grid.
std::vector< bool > periodic
 Whether some colvars are periodic.
std::vector< cvm::realwidths
 Widths of the colvars in this grid.
bool has_parent_data
 True if this is a count grid related to another grid of data.

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.

Detailed Description

template<class T>
class colvar_grid< T >

Grid of values of a function of several collective variables.

Parameters:
T The data type
Only scalar colvars supported so far

Definition at line 19 of file colvargrid.h.


Constructor & Destructor Documentation

template<class T>
colvar_grid< T >::colvar_grid  )  [inline]
 

Default constructor.

Definition at line 145 of file colvargrid.h.

00146   {
00147     save_delimiters = false;
00148     nd = nt = 0;
00149   }

template<class T>
virtual colvar_grid< T >::~colvar_grid  )  [inline, virtual]
 

Destructor.

Definition at line 152 of file colvargrid.h.

00153   {}

template<class T>
colvar_grid< T >::colvar_grid colvar_grid< T > const &  g  )  [inline]
 

"Almost copy-constructor": only copies configuration parameters from another grid, but doesn't reallocate stuff; create() must be called after that;

Definition at line 158 of file colvargrid.h.

00159     : nd (g.nd),
00160       mult (g.mult),
00161       cv (g.cv),
00162       lower_boundaries (g.lower_boundaries),
00163       upper_boundaries (g.upper_boundaries),
00164       periodic (g.periodic),
00165       widths (g.widths),
00166       data()
00167   {
00168     save_delimiters = false;
00169   }

template<class T>
colvar_grid< T >::colvar_grid std::vector< int > const &  nx_i,
T const &  t = T(),
size_t const &  mult_i = 1
[inline]
 

Constructor from explicit grid sizes.

Parameters:
nx_i Number of grid points along each dimension
t Initial value for the function at each point (optional)
mult_i Multiplicity of each value

Definition at line 175 of file colvargrid.h.

00178   {
00179     save_delimiters = false;
00180     this->create (nx_i, t, mult_i);
00181   }

template<class T>
colvar_grid< T >::colvar_grid std::vector< colvar * > const &  colvars,
T const &  t = T(),
size_t const &  mult_i = 1,
bool  margin = false
[inline]
 

Constructor from a vector of colvars.

Definition at line 184 of file colvargrid.h.

00188     : cv (colvars)
00189   {
00190     save_delimiters = false;
00191 
00192     std::vector<int> nx_i;
00193 
00194     if (cvm::debug())
00195       cvm::log ("Allocating a grid for "+cvm::to_str (colvars.size())+
00196                 " collective variables.\n");
00197 
00198     for (size_t i =  0; i < cv.size(); i++) {
00199 
00200       if (cv[i]->type() != colvarvalue::type_scalar) {
00201         cvm::fatal_error ("Colvar grids can only be automatically "
00202                           "constructed for scalar variables.  "
00203                           "ABF and histogram can not be used; metadynamics "
00204                           "can be used with useGrids disabled.\n");
00205       }
00206 
00207       if (cv[i]->width <= 0.0) {
00208         cvm::fatal_error ("Tried to initialize a grid on a"
00209                           "variable with negative or zero width.\n");
00210       }
00211 
00212       if (!cv[i]->tasks[colvar::task_lower_boundary] || !cv[i]->tasks[colvar::task_upper_boundary]) {
00213         cvm::fatal_error ("Tried to initialize a grid on a"
00214                           "variable with undefined boundaries.\n");
00215       }
00216 
00217       widths.push_back (cv[i]->width);
00218       periodic.push_back (cv[i]->periodic_boundaries());
00219 
00220       if (margin) {
00221         if (periodic[i]) {
00222           // Shift the grid by half the bin width (values at edges instead of center of bins)
00223           lower_boundaries.push_back (cv[i]->lower_boundary.real_value - 0.5 * widths[i]);
00224           upper_boundaries.push_back (cv[i]->upper_boundary.real_value - 0.5 * widths[i]);
00225         } else {
00226           // Make this grid larger by one bin width
00227           lower_boundaries.push_back (cv[i]->lower_boundary.real_value - 0.5 * widths[i]);
00228           upper_boundaries.push_back (cv[i]->upper_boundary.real_value + 0.5 * widths[i]);
00229         }
00230       } else {
00231         lower_boundaries.push_back (cv[i]->lower_boundary);
00232         upper_boundaries.push_back (cv[i]->upper_boundary);
00233       }
00234 
00235 
00236       { 
00237         cvm::real nbins = ( upper_boundaries[i].real_value -
00238                             lower_boundaries[i].real_value ) / widths[i];
00239         //int nbins_round = ::round (nbins);
00240         int nbins_round = (int)(nbins+0.5);
00241         
00242         if (::fabs (nbins - cvm::real (nbins_round)) > 1.0E-10) {
00243           cvm::log ("Warning: grid interval ("+
00244                     cvm::to_str (lower_boundaries[i], cvm::cv_width, cvm::cv_prec)+" - "+
00245                     cvm::to_str (upper_boundaries[i], cvm::cv_width, cvm::cv_prec)+
00246                     ") is not commensurate to its bin width ("+
00247                     cvm::to_str (widths[i], cvm::cv_width, cvm::cv_prec)+").\n");
00248           upper_boundaries[i].real_value = lower_boundaries[i].real_value +
00249             (nbins_round * widths[i]);
00250         }
00251 
00252         if (cvm::debug())
00253           cvm::log ("Number of points is "+cvm::to_str ((int) nbins_round)+
00254                     " for the colvar no. "+cvm::to_str (i+1)+".\n");
00255         
00256         nx_i.push_back (nbins_round);
00257       }
00258 
00259     }
00260 
00261     create (nx_i, t, mult_i);
00262   }


Member Function Documentation

template<class T>
void colvar_grid< T >::add_constant T const &  t  )  [inline]
 

Add a constant to all elements (fast loop).

Definition at line 315 of file colvargrid.h.

Referenced by colvarbias_meta::write_restart().

00316   {
00317     for (size_t i = 0; i < nt; i++) 
00318       data[i] += t;
00319   }

template<class T>
size_t colvar_grid< T >::address std::vector< int > const &  ix  )  const [inline, protected]
 

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   }

template<class T>
cvm::real colvar_grid< T >::bin_distance_from_boundaries std::vector< colvarvalue > const &  values  )  [inline]
 

Get the minimum distance (in number of bins) from the boundaries; give a negative number if the point is off-grid.

Definition at line 353 of file colvargrid.h.

Referenced by colvarbias_meta::create_hill(), colvarbias_meta::read_hill(), and colvarbias_meta::recount_hills_off_grid().

00354   {
00355     cvm::real minimum = 1.0E+16;
00356     for (size_t i = 0; i < nd; i++) {
00357 
00358       cvm::real dl = ::sqrt (cv[i]->dist2 (values[i], lower_boundaries[i])) / widths[i];
00359       cvm::real du = ::sqrt (cv[i]->dist2 (values[i], upper_boundaries[i])) / widths[i];
00360 
00361       if (! periodic[i]) {
00362         if (values[i].real_value < lower_boundaries[i])
00363           dl *= -1.0;
00364         if (values[i].real_value > upper_boundaries[i])
00365           du *= -1.0;
00366       }
00367 
00368       if (dl < minimum) 
00369         minimum = dl;
00370       if (du < minimum) 
00371         minimum = du;
00372     }
00373 
00374     return minimum;
00375   }

template<class T>
colvarvalue colvar_grid< T >::bin_to_value_scalar int const &  i_bin,
int const   i
const [inline]
 

Use the two boundaries and the width to report the central value corresponding to a bin index.

Definition at line 291 of file colvargrid.h.

Referenced by colvarbias_meta::project_hills(), and colvar_grid< size_t >::write_multicol().

00292   {
00293     return lower_boundaries[i].real_value + widths[i] * (0.5 + i_bin);
00294   }

template<class T>
void colvar_grid< T >::check_consistency  )  [inline]
 

Check that the grid information inside (boundaries, widths, ...) is consistent with the current setting of the colvars.

Definition at line 571 of file colvargrid.h.

00572   {
00573     for (size_t i = 0; i < nd; i++) {
00574       if ( (::sqrt (cv[i]->dist2 (cv[i]->lower_boundary,
00575                                   lower_boundaries[i])) > 1.0E-10) || 
00576            (::sqrt (cv[i]->dist2 (cv[i]->upper_boundary,
00577                                   upper_boundaries[i])) > 1.0E-10) || 
00578            (::sqrt (cv[i]->dist2 (cv[i]->width,
00579                                   widths[i])) > 1.0E-10) ) {
00580         cvm::fatal_error ("Error: restart information for a grid is "
00581                           "inconsistent with that of its colvars.\n");
00582       }
00583     }
00584   }

template<class T>
void colvar_grid< T >::create  )  [inline]
 

Allocate data (allow initialization also after construction).

Definition at line 138 of file colvargrid.h.

Referenced by colvar_grid< size_t >::colvar_grid(), colvar_grid< size_t >::create(), and colvar_grid< size_t >::parse_params().

00139   {
00140     create (this->nx, T(), this->mult);
00141   }

template<class T>
void colvar_grid< T >::create std::vector< int > const &  nx_i,
T const &  t = T(),
size_t const &  mult_i = 1
[inline]
 

Allocate data (allow initialization also after construction).

Definition at line 115 of file colvargrid.h.

Referenced by colvarbias_meta::update().

00118   {
00119     mult = mult_i;
00120     nd = nx_i.size();
00121     nxc.resize (nd);
00122     nx = nx_i;
00123 
00124     nt = mult;
00125     for (int i = nd-1; i >= 0; i--) {
00126       if (nx_i[i] <= 0)
00127         cvm::fatal_error ("Error: providing an invalid number of points, "+
00128                           cvm::to_str (nx_i[i])+".\n");
00129       nxc[i] = nt;
00130       nt *= nx[i];
00131     }
00132 
00133     data.reserve (nt);
00134     data.assign (nt, t);
00135   }

template<class T>
std::vector<int> const colvar_grid< T >::get_colvars_index  )  const [inline]
 

Get the bin indices corresponding to the current values of the colvars.

Definition at line 342 of file colvargrid.h.

Referenced by colvar_grid< size_t >::get_colvars_index().

00343   {
00344     std::vector<colvarvalue> current_values (nd);
00345     for (size_t i = 0; i < nd; i++) {
00346       current_values[i] = cv[i]->value();
00347     }
00348     return get_colvars_index (current_values);
00349   }

template<class T>
std::vector<int> const colvar_grid< T >::get_colvars_index std::vector< colvarvalue > const &  values  )  const [inline]
 

Get the bin indices corresponding to the provided values of the colvars.

Definition at line 331 of file colvargrid.h.

Referenced by colvarbias_meta::update().

00332   {
00333     std::vector<int> index = new_index();
00334     for (size_t i = 0; i < nd; i++) {
00335       index[i] = cv[i]->value_to_bin_scalar (values[i], lower_boundaries[i], widths[i]);
00336     }
00337     return index;
00338   }

template<class T>
void colvar_grid< T >::incr std::vector< int > &  ix  )  const [inline]
 

Increment the index, in a way that will make it loop over the whole nd-dimensional array.

Definition at line 470 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().

00471   {
00472     for (int i = ix.size()-1; i >= 0; i--) {
00473 
00474       ix[i]++;
00475 
00476       if (ix[i] >= nx[i]) {
00477 
00478         if (i > 0) {
00479           ix[i] = 0;
00480           continue;
00481         } else {
00482           // this is the last iteration, a non-valid index is being
00483           // set for the outer index, which will be caught by
00484           // index_ok()
00485           ix[0] = nx[0];
00486           return;
00487         }
00488       } else {
00489         return;
00490       }
00491     }
00492   }

template<class T>
bool colvar_grid< T >::index_ok std::vector< int > const &  ix  )  const [inline]
 

Check that the index is within range in each of the dimensions.

Definition at line 459 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().

00460   {
00461     for (size_t i = 0; i < nd; i++) {
00462       if ( (ix[i] < 0) || (ix[i] >= int (nx[i])) )
00463         return false;
00464     }
00465     return true;
00466   }

template<class T>
void colvar_grid< T >::map_grid colvar_grid< T > const &  other_grid  )  [inline]
 

Add data from another grid of the same type.

Parameters:
grid Reference to the other grid
other_grid_boundaries Lower boundaries in the other grid (by default same as this)
other_grid_widths in the other grid (by default same as this)
Note: this function maps other_grid inside this one regardless of whether it fits or not.

Definition at line 385 of file colvargrid.h.

Referenced by colvarbias_meta::read_restart(), and colvarbias_meta::update().

00386   {
00387     if (other_grid.multiplicity() != this->multiplicity())
00388       cvm::fatal_error ("Error: trying to merge two grids with values of "
00389                         "different multiplicity.\n");
00390 
00391     std::vector<colvarvalue> const &gb  = this->lower_boundaries;
00392     std::vector<cvm::real> const &gw    = this->widths;
00393     std::vector<colvarvalue> const &ogb = other_grid.lower_boundaries;
00394     std::vector<cvm::real> const &ogw   = other_grid.widths;
00395 
00396     std::vector<int> ix = this->new_index();
00397     std::vector<int> oix = other_grid.new_index();
00398 
00399     if (cvm::debug())
00400       cvm::log ("Remapping grid...\n");
00401     for ( ; this->index_ok (ix); this->incr (ix)) {
00402 
00403       for (size_t i = 0; i < nd; i++) {
00404         oix[i] =
00405           cv[i]->value_to_bin_scalar (cv[i]->bin_to_value_scalar (ix[i], gb[i], gw[i]),
00406                                       ogb[i],
00407                                       ogw[i]);
00408       }
00409 
00410       if (! other_grid.index_ok (oix)) {
00411         continue;
00412       }
00413 
00414       for (size_t im = 0; im < mult; im++) {
00415         this->set_value (ix, other_grid.value (oix, im), im);
00416       }
00417     }
00418     if (cvm::debug())
00419       cvm::log ("Remapping done.\n");
00420   }

template<class T>
size_t colvar_grid< T >::multiplicity  )  const [inline]
 

Return the multiplicity of the type used.

Definition at line 109 of file colvargrid.h.

Referenced by colvar_grid< size_t >::map_grid().

00110   {
00111     return mult;
00112   }

template<class T>
void colvar_grid< T >::multiply_constant cvm::real const &  a  )  [inline]
 

Multiply all elements by a scalar constant (fast loop).

Definition at line 322 of file colvargrid.h.

Referenced by colvarbias_meta::write_restart().

00323   {
00324     for (size_t i = 0; i < nt; i++) 
00325       data[i] *= a;
00326   }

template<class T>
std::vector<int> const colvar_grid< T >::new_index  )  const [inline]
 

Get the index corresponding to the "first" bin, to be used as the initial value for an index in looping.

Definition at line 452 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().

00453   {
00454     return std::vector<int> (nd, 0);
00455   }

template<class T>
size_t colvar_grid< T >::number_of_colvars  )  const [inline]
 

Return the number of colvars.

Definition at line 80 of file colvargrid.h.

00081   {
00082     return nd;
00083   }

template<class T>
size_t colvar_grid< T >::number_of_points int const   icv = -1  )  const [inline]
 

Return the number of points in the i-th direction, if provided, or the total number

Definition at line 87 of file colvargrid.h.

Referenced by colvar_grid< size_t >::read_raw().

00088   {
00089     if (icv < 0) {
00090       return nt;
00091     } else {
00092       return nx[icv];
00093     }
00094   }

template<class T>
bool colvar_grid< T >::parse_params std::string const &  conf  )  [inline]
 

Definition at line 524 of file colvargrid.h.

00525   {
00526     std::vector<int> old_nx = nx;
00527     std::vector<colvarvalue> old_lb = lower_boundaries;
00528 
00529     {
00530       size_t nd_in = 0;
00531       colvarparse::get_keyval (conf, "n_colvars", nd_in, nd, colvarparse::parse_silent);
00532       if (nd_in != nd)
00533         cvm::fatal_error ("Error: trying to read data for a grid "
00534                           "that contains a different number of colvars ("+
00535                           cvm::to_str (nd_in)+") than the grid defined "
00536                           "in the configuration file ("+cvm::to_str (nd)+
00537                           ").\n");
00538     }
00539 
00540     colvarparse::get_keyval (conf, "lower_boundaries",
00541                              lower_boundaries, lower_boundaries, colvarparse::parse_silent);
00542 
00543     colvarparse::get_keyval (conf, "upper_boundaries",
00544                              upper_boundaries, upper_boundaries, colvarparse::parse_silent);
00545 
00546     colvarparse::get_keyval (conf, "widths", widths, widths, colvarparse::parse_silent);
00547 
00548     colvarparse::get_keyval (conf, "sizes", nx, nx, colvarparse::parse_silent);
00549 
00550     bool new_params = false;
00551     for (size_t i = 0; i < nd; i++) {
00552       if ( (old_nx[i] != nx[i]) ||
00553            (::sqrt (cv[i]->dist2 (old_lb[i],
00554                                   lower_boundaries[i])) > 1.0E-10) ) {
00555         new_params = true;
00556       }
00557     }
00558 
00559     // reallocate the array in case the grid params have just changed
00560     if (new_params) {
00561       data.resize (0);
00562       this->create (nx, T(), mult);
00563     }
00564 
00565     return true;
00566   }

template<class T>
void colvar_grid< T >::read_multicol std::istream &  is,
bool  add = false
[inline]
 

Read a grid written by colvar_grid::write_multicol() Adding data if add is true, replacing if false.

Definition at line 684 of file colvargrid.h.

00685   {
00686     // Data in the header: nColvars, then for each
00687     // xiMin, dXi, nPoints, periodic
00688 
00689     std::string   hash;
00690     cvm::real     lower, width, x;
00691     size_t        n, periodic;
00692     bool          remap;
00693     std::vector<T>        new_value;
00694     std::vector<int>      nx_read;
00695     std::vector<int>      bin;
00696 
00697     if ( cv.size() != nd ) {
00698       cvm::fatal_error ("Cannot read grid file: missing reference to colvars.");
00699     }
00700 
00701     if ( !(is >> hash) || (hash != "#") ) {
00702       cvm::fatal_error ("Error reading grid at position "+
00703                         cvm::to_str (is.tellg())+" in stream (read \"" + hash + "\")\n");
00704     }
00705 
00706     is >> n;
00707     if ( n != nd ) {
00708       cvm::fatal_error ("Error reading grid: wrong number of collective variables.\n");
00709     }
00710 
00711     nx_read.resize (n);
00712     bin.resize (n);
00713     new_value.resize (mult);
00714 
00715     if (this->has_parent_data && add) {
00716       new_data.resize (data.size());
00717     }
00718 
00719     remap = false;
00720     for (size_t i = 0; i < nd; i++ ) {
00721       if ( !(is >> hash) || (hash != "#") ) {
00722         cvm::fatal_error ("Error reading grid at position "+
00723                           cvm::to_str (is.tellg())+" in stream (read \"" + hash + "\")\n");
00724       }
00725 
00726       is >> lower >> width >> nx_read[i] >> periodic;
00727 
00728 
00729       if ( (::fabs (lower - lower_boundaries[i].real_value) > 1.0e-10) ||
00730            (::fabs (width - widths[i] ) > 1.0e-10) ||
00731            (nx_read[i] != nx[i]) ) {
00732         cvm::log ("Warning: reading from different grid definition (colvar "
00733                   + cvm::to_str (i+1) + "); remapping data on new grid.\n");
00734         remap = true;
00735       }
00736     }
00737 
00738     if ( remap ) {
00739       // re-grid data
00740       while (is.good()) {
00741         bool end_of_file = false;
00742 
00743         for (size_t i = 0; i < nd; i++ ) {
00744           if ( !(is >> x) ) end_of_file = true;
00745           bin[i] = value_to_bin_scalar (x, i);
00746         }
00747         if (end_of_file) break; 
00748 
00749         for (size_t imult = 0; imult < mult; imult++) {
00750           is >> new_value[imult];
00751         }
00752 
00753         if ( index_ok(bin) ) {
00754           for (size_t imult = 0; imult < mult; imult++) {
00755             value_input (bin, new_value[imult], imult, add);
00756           }
00757         }
00758       }
00759     } else {
00760       // do not re-grid the data but assume the same grid is used
00761       for (std::vector<int> ix = new_index(); index_ok (ix); incr (ix) ) {
00762         for (size_t i = 0; i < nd; i++ ) {
00763           is >> x;
00764         }
00765         for (size_t imult = 0; imult < mult; imult++) {
00766           is >> new_value[imult];
00767           value_input (ix, new_value[imult], imult, add);
00768         }
00769       }
00770     }
00771     return;
00772   }

template<class T>
std::istream& colvar_grid< T >::read_raw std::istream &  is  )  [inline]
 

Read data written by colvar_grid::write_raw().

Definition at line 615 of file colvargrid.h.

00616   {
00617     std::vector<int> ix = new_index();
00618     size_t count = 0;
00619     for ( ; index_ok (ix); incr (ix)) {
00620       for (size_t imult = 0; imult < mult; imult++) {
00621         T new_value;
00622         if (is >> new_value) {
00623           value_input (ix, new_value, imult);
00624           count++;
00625         } else {
00626           cvm::fatal_error ("Error: malformed grid data in input.\n");
00627         }
00628       }
00629     }
00630 
00631     if (count != number_of_points()) {
00632       cvm::log ("Error: number of grid points read from file ("+
00633                 cvm::to_str (count)+") is wrong (should be "+
00634                 cvm::to_str (number_of_points())+").\n");
00635       cvm::fatal_error ("Possible explanations: grid parameters in the configuration (lowerBoundary, upperBoundary, width) are different from those in the file, or the file is corrupt.\n");
00636     }
00637 
00638     return is;
00639   }

template<class T>
void colvar_grid< T >::set_sizes std::vector< int > const &  new_sizes  )  [inline]
 

Set the sizes in each direction.

Definition at line 103 of file colvargrid.h.

00104   {
00105     nx = new_sizes;
00106   }

template<class T>
void colvar_grid< T >::set_value std::vector< int > const &  ix,
T const &  t,
size_t const &  imult = 0
[inline]
 

Set the value at the point with index ix.

Definition at line 297 of file colvargrid.h.

Referenced by colvar_grid< size_t >::map_grid().

00300   {
00301     data[this->address (ix)+imult] = t;
00302   }

template<class T>
std::vector<int> const& colvar_grid< T >::sizes  )  const [inline]
 

Get the sizes in each direction.

Definition at line 97 of file colvargrid.h.

Referenced by colvarbias_meta::update().

00098   {
00099     return nx;
00100   }

template<class T>
T const& colvar_grid< T >::value std::vector< int > const &  ix,
size_t const &  imult = 0
const [inline]
 

Get the binned value indexed by ix, or the first of them if the multiplicity is larger than 1.

Definition at line 307 of file colvargrid.h.

Referenced by colvar_grid_gradient::average(), colvarbias_meta::calc_hills(), 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().

00309   {
00310     return data[this->address (ix) + imult];
00311   }

template<class T>
virtual void colvar_grid< T >::value_input std::vector< int > const &  ix,
T const &  t,
size_t const &  imult = 0,
bool  add = false
[inline, virtual]
 

Get the value from a formatted output and transform it into the internal representation (it may have been rescaled or manipulated).

Reimplemented in colvar_grid_count, colvar_grid_scalar, and colvar_grid_gradient.

Definition at line 433 of file colvargrid.h.

Referenced by colvar_grid< size_t >::read_multicol(), and colvar_grid< size_t >::read_raw().

00437   {
00438     if ( add )
00439       data[address (ix) + imult] += t;
00440     else
00441       data[address (ix) + imult] = t;
00442   }

template<class T>
virtual T colvar_grid< T >::value_output std::vector< int > const &  ix,
size_t const &  imult = 0
[inline, virtual]
 

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 424 of file colvargrid.h.

Referenced by colvar_grid< size_t >::write_multicol(), and colvar_grid< size_t >::write_raw().

00426   {
00427     return value (ix, imult);
00428   }

template<class T>
int colvar_grid< T >::value_to_bin_scalar colvarvalue const &  value,
const int  i
const [inline]
 

Use the lower boundary and the width to report which bin the provided value is in.

Definition at line 284 of file colvargrid.h.

Referenced by colvar_grid< size_t >::map_grid(), and colvar_grid< size_t >::read_multicol().

00285   {
00286     return (int) ::floor ( (value.real_value - lower_boundaries[i].real_value) / widths[i] );
00287   }

template<class T>
void colvar_grid< T >::wrap std::vector< int > &  ix  )  const [inline]
 

Wrap an index vector around periodic boundary conditions also checks validity of non-periodic indices

Definition at line 267 of file colvargrid.h.

00268   {
00269     for (size_t i = 0; i < nd; i++) {
00270       if (periodic[i]) {
00271         ix[i] = (ix[i] + nx[i]) % nx[i]; //to ensure non-negative result
00272       } else {
00273         if (ix[i] < 0 || ix[i] >= nx[i])
00274           cvm::fatal_error ("Trying to wrap illegal index vector (non-PBC): "
00275                             + cvm::to_str (ix));
00276       }
00277     }
00278   }

template<class T>
void colvar_grid< T >::write_multicol std::ostream &  os  )  [inline]
 

Write the grid in a format which is both human readable and suitable for visualization e.g. with gnuplot.

Definition at line 643 of file colvargrid.h.

Referenced by colvarbias_histogram::update(), and colvarbias_meta::write_restart().

00644   {
00645     std::streamsize const w = os.width();
00646     std::streamsize const p = os.precision();
00647 
00648     // Data in the header: nColvars, then for each
00649     // xiMin, dXi, nPoints, periodic
00650 
00651     os << std::setw (2) << "# " << nd << "\n";
00652     for (size_t i = 0; i < nd; i++) {
00653       os << "# "
00654          << std::setw (10) << lower_boundaries[i]
00655          << std::setw (10) << widths[i]
00656          << std::setw (10) << nx[i] << "  "
00657          << periodic[i] << "\n";
00658     }
00659 
00660     for (std::vector<int> ix = new_index(); index_ok (ix); incr (ix) ) {
00661 
00662       if (ix.back() == 0) {
00663         // if the last index is 0, add a new line to mark the new record
00664         os << "\n";
00665       }
00666 
00667       for (size_t i = 0; i < nd; i++) {
00668         os << " "
00669            << std::setw (w) << std::setprecision (p)
00670            << bin_to_value_scalar (ix[i], i);
00671       }
00672       os << " ";
00673       for (size_t imult = 0; imult < mult; imult++) {
00674         os << " "
00675            << std::setw (w) << std::setprecision (p)
00676            << value_output (ix, imult);
00677       }
00678       os << "\n";
00679     }
00680   }

template<class T>
std::ostream& colvar_grid< T >::write_params std::ostream &  os  )  [inline]
 

Write the grid parameters (number of colvars, boundaries, width and number of points).

Definition at line 495 of file colvargrid.h.

00496   {
00497     os << "grid_parameters {\n  n_colvars " << nd << "\n";
00498       
00499     os << "  lower_boundaries ";
00500     for (size_t i = 0; i < nd; i++) 
00501       os << " " << lower_boundaries[i];
00502     os << "\n";
00503 
00504     os << "  upper_boundaries ";
00505     for (size_t i = 0; i < nd; i++) 
00506       os << " " << upper_boundaries[i];
00507     os << "\n";
00508 
00509     os << "  widths ";
00510     for (size_t i = 0; i < nd; i++) 
00511       os << " " << widths[i];
00512     os << "\n";
00513 
00514     os << "  sizes ";
00515     for (size_t i = 0; i < nd; i++) 
00516       os << " " << nx[i];
00517     os << "\n";
00518 
00519     os << "}\n";
00520     return os;
00521   } 

template<class T>
std::ostream& colvar_grid< T >::write_raw std::ostream &  os,
size_t const   buf_size = 3
[inline]
 

Write the grid data without labels, as they are represented in memory.

Parameters:
buf_size Number of values per line

Definition at line 590 of file colvargrid.h.

00592   {
00593     std::streamsize const w = os.width();
00594     std::streamsize const p = os.precision();
00595 
00596     std::vector<int> ix = new_index();
00597     size_t count = 0;
00598     for ( ; index_ok (ix); incr (ix)) {
00599       for (size_t imult = 0; imult < mult; imult++) {
00600         os << " "
00601            << std::setw (w) << std::setprecision (p)
00602            << value_output (ix, imult);
00603         if (((++count) % buf_size) == 0)
00604           os << "\n";
00605       }
00606     }
00607     // write a final newline only if buffer is not empty
00608     if ((count % buf_size) != 0)
00609       os << "\n";
00610 
00611     return os;
00612   }


Member Data Documentation

template<class T>
std::vector<colvar *> colvar_grid< T >::cv [protected]
 

Colvars collected in this grid.

Definition at line 46 of file colvargrid.h.

template<class T>
std::vector<T> colvar_grid< T >::data [protected]
 

Low-level array of values.

Definition at line 40 of file colvargrid.h.

template<class T>
bool colvar_grid< T >::has_parent_data
 

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().

template<class T>
std::vector<colvarvalue> colvar_grid< T >::lower_boundaries
 

Lower boundaries of the colvars in this grid.

Definition at line 65 of file colvargrid.h.

Referenced by colvar_grid< size_t >::map_grid(), and colvarbias_meta::update().

template<class T>
size_t colvar_grid< T >::mult [protected]
 

Multiplicity of each datum (allow the binning of non-scalar types).

Definition at line 34 of file colvargrid.h.

template<class T>
size_t colvar_grid< T >::nd [protected]
 

Number of dimensions.

Definition at line 24 of file colvargrid.h.

template<class T>
std::vector<size_t> colvar_grid< T >::new_data [protected]
 

Newly read data (used for count grids, when adding several grids read from disk).

Definition at line 43 of file colvargrid.h.

template<class T>
size_t colvar_grid< T >::nt [protected]
 

Total number of grid points.

Definition at line 37 of file colvargrid.h.

template<class T>
std::vector<int> colvar_grid< T >::nx [protected]
 

Number of points along each dimension.

Definition at line 27 of file colvargrid.h.

template<class T>
std::vector<int> colvar_grid< T >::nxc [protected]
 

Cumulative number of points along each dimension.

Definition at line 30 of file colvargrid.h.

template<class T>
std::vector<bool> colvar_grid< T >::periodic
 

Whether some colvars are periodic.

Definition at line 71 of file colvargrid.h.

template<class T>
std::vector<colvarvalue> colvar_grid< T >::upper_boundaries
 

Upper boundaries of the colvars in this grid.

Definition at line 68 of file colvargrid.h.

Referenced by colvarbias_meta::update().

template<class T>
std::vector<cvm::real> colvar_grid< T >::widths
 

Widths of the colvars in this grid.

Definition at line 74 of file colvargrid.h.

Referenced by colvarbias_meta::calc_hills(), and colvar_grid< size_t >::map_grid().


The documentation for this class was generated from the following file:
Generated on Mon Nov 23 04:59:32 2009 for NAMD by  doxygen 1.3.9.1