2016-06-18 07:07:51 +08:00
// -*- c++ -*-
2012-05-24 00:20:04 +08:00
# ifndef COLVARGRID_H
# define COLVARGRID_H
# include <iostream>
# include <iomanip>
# include <cmath>
# include "colvar.h"
# include "colvarmodule.h"
# include "colvarvalue.h"
# include "colvarparse.h"
/// \brief Grid of values of a function of several collective
/// variables \param T The data type
///
2015-02-20 02:20:52 +08:00
/// Only scalar colvars supported so far: vector colvars are treated as arrays
2012-05-24 00:20:04 +08:00
template < class T > class colvar_grid : public colvarparse {
protected :
/// Number of dimensions
size_t nd ;
/// Number of points along each dimension
std : : vector < int > nx ;
/// Cumulative number of points along each dimension
std : : vector < int > nxc ;
/// \brief Multiplicity of each datum (allow the binning of
2014-12-02 10:09:53 +08:00
/// non-scalar types such as atomic gradients)
2012-05-24 00:20:04 +08:00
size_t mult ;
/// Total number of grid points
size_t nt ;
/// Low-level array of values
std : : vector < T > data ;
/// Newly read data (used for count grids, when adding several grids read from disk)
std : : vector < size_t > new_data ;
/// Colvars collected in this grid
std : : vector < colvar * > cv ;
2012-06-30 01:52:31 +08:00
/// Do we request actual value (for extended-system colvars)?
std : : vector < bool > actual_value ;
2012-05-24 00:20:04 +08:00
/// Get the low-level index corresponding to an index
2014-12-02 10:09:53 +08:00
inline size_t address ( std : : vector < int > const & ix ) const
2012-05-24 00:20:04 +08:00
{
size_t addr = 0 ;
for ( size_t i = 0 ; i < nd ; i + + ) {
addr + = ix [ i ] * nxc [ i ] ;
if ( cvm : : debug ( ) ) {
2014-10-08 04:30:25 +08:00
if ( ix [ i ] > = nx [ i ] ) {
2015-02-20 02:20:52 +08:00
cvm : : error ( " Error: exceeding bounds in colvar_grid. \n " , BUG_ERROR ) ;
2014-10-08 04:30:25 +08:00
return 0 ;
}
2012-05-24 00:20:04 +08:00
}
}
return addr ;
}
public :
/// Lower boundaries of the colvars in this grid
std : : vector < colvarvalue > lower_boundaries ;
/// Upper boundaries of the colvars in this grid
std : : vector < colvarvalue > upper_boundaries ;
/// Whether some colvars are periodic
std : : vector < bool > periodic ;
2012-06-30 01:52:31 +08:00
/// Whether some colvars have hard lower boundaries
std : : vector < bool > hard_lower_boundaries ;
/// Whether some colvars have hard upper boundaries
std : : vector < bool > hard_upper_boundaries ;
2012-05-24 00:20:04 +08:00
/// Widths of the colvars in this grid
std : : vector < cvm : : real > widths ;
/// True if this is a count grid related to another grid of data
bool has_parent_data ;
/// Whether this grid has been filled with data or is still empty
bool has_data ;
/// Return the number of colvars
inline size_t number_of_colvars ( ) const
{
return nd ;
}
/// Return the number of points in the i-th direction, if provided, or
/// the total number
2014-12-02 10:09:53 +08:00
inline size_t number_of_points ( int const icv = - 1 ) const
2012-05-24 00:20:04 +08:00
{
if ( icv < 0 ) {
return nt ;
} else {
return nx [ icv ] ;
}
}
/// Get the sizes in each direction
inline std : : vector < int > const & sizes ( ) const
{
return nx ;
}
/// Set the sizes in each direction
2014-12-02 10:09:53 +08:00
inline void set_sizes ( std : : vector < int > const & new_sizes )
2012-05-24 00:20:04 +08:00
{
nx = new_sizes ;
}
/// Return the multiplicity of the type used
inline size_t multiplicity ( ) const
{
return mult ;
}
2016-08-04 03:32:03 +08:00
/// \brief Request grid to use actual values of extended coords
inline void request_actual_value ( bool b = true )
{
size_t i ;
for ( i = 0 ; i < actual_value . size ( ) ; i + + ) {
actual_value [ i ] = b ;
}
}
2015-02-20 02:20:52 +08:00
/// \brief Allocate data
int setup ( std : : vector < int > const & nx_i ,
T const & t = T ( ) ,
size_t const & mult_i = 1 )
2012-05-24 00:20:04 +08:00
{
2015-02-20 02:20:52 +08:00
if ( cvm : : debug ( ) ) {
cvm : : log ( " Allocating grid: multiplicity = " + cvm : : to_str ( mult_i ) +
" , dimensions = " + cvm : : to_str ( nx_i ) + " . \n " ) ;
}
2012-05-24 00:20:04 +08:00
mult = mult_i ;
2015-02-20 02:20:52 +08:00
data . clear ( ) ;
2012-05-24 00:20:04 +08:00
nx = nx_i ;
2015-02-20 02:20:52 +08:00
nd = nx . size ( ) ;
nxc . resize ( nd ) ;
2012-05-24 00:20:04 +08:00
2015-02-20 02:20:52 +08:00
// setup dimensions
2012-05-24 00:20:04 +08:00
nt = mult ;
for ( int i = nd - 1 ; i > = 0 ; i - - ) {
2015-02-20 02:20:52 +08:00
if ( nx [ i ] < = 0 ) {
cvm : : error ( " Error: providing an invalid number of grid points, " +
cvm : : to_str ( nx [ i ] ) + " . \n " , BUG_ERROR ) ;
return COLVARS_ERROR ;
2014-10-08 04:30:25 +08:00
}
2012-05-24 00:20:04 +08:00
nxc [ i ] = nt ;
nt * = nx [ i ] ;
}
2015-02-20 02:20:52 +08:00
if ( cvm : : debug ( ) ) {
cvm : : log ( " Total number of grid elements = " + cvm : : to_str ( nt ) + " . \n " ) ;
}
2014-12-02 10:09:53 +08:00
data . reserve ( nt ) ;
data . assign ( nt , t ) ;
2015-02-20 02:20:52 +08:00
return COLVARS_OK ;
2012-05-24 00:20:04 +08:00
}
/// \brief Allocate data (allow initialization also after construction)
2015-02-20 02:20:52 +08:00
int setup ( )
2012-05-24 00:20:04 +08:00
{
2015-02-20 02:20:52 +08:00
return setup ( this - > nx , T ( ) , this - > mult ) ;
2012-05-24 00:20:04 +08:00
}
/// \brief Reset data (in case the grid is being reused)
2014-12-02 10:09:53 +08:00
void reset ( T const & t = T ( ) )
2012-05-24 00:20:04 +08:00
{
2014-12-02 10:09:53 +08:00
data . assign ( nt , t ) ;
2012-05-24 00:20:04 +08:00
}
/// Default constructor
2014-12-02 10:09:53 +08:00
colvar_grid ( ) : has_data ( false )
2012-05-24 00:20:04 +08:00
{
save_delimiters = false ;
nd = nt = 0 ;
2015-02-20 02:20:52 +08:00
mult = 1 ;
this - > setup ( ) ;
2012-05-24 00:20:04 +08:00
}
/// Destructor
virtual ~ colvar_grid ( )
{ }
/// \brief "Almost copy-constructor": only copies configuration
/// parameters from another grid, but doesn't reallocate stuff;
2015-02-20 02:20:52 +08:00
/// setup() must be called after that;
2014-12-02 10:09:53 +08:00
colvar_grid ( colvar_grid < T > const & g ) : nd ( g . nd ) ,
2015-02-20 02:20:52 +08:00
nx ( g . nx ) ,
mult ( g . mult ) ,
data ( ) ,
cv ( g . cv ) ,
actual_value ( g . actual_value ) ,
lower_boundaries ( g . lower_boundaries ) ,
upper_boundaries ( g . upper_boundaries ) ,
periodic ( g . periodic ) ,
hard_lower_boundaries ( g . hard_lower_boundaries ) ,
hard_upper_boundaries ( g . hard_upper_boundaries ) ,
widths ( g . widths ) ,
has_data ( false )
2012-05-24 00:20:04 +08:00
{
save_delimiters = false ;
}
/// \brief Constructor from explicit grid sizes \param nx_i Number
/// of grid points along each dimension \param t Initial value for
/// the function at each point (optional) \param mult_i Multiplicity
/// of each value
2014-12-02 10:09:53 +08:00
colvar_grid ( std : : vector < int > const & nx_i ,
2015-02-20 02:20:52 +08:00
T const & t = T ( ) ,
size_t mult_i = 1 )
: has_data ( false )
2012-05-24 00:20:04 +08:00
{
save_delimiters = false ;
2015-02-20 02:20:52 +08:00
this - > setup ( nx_i , t , mult_i ) ;
2012-05-24 00:20:04 +08:00
}
/// \brief Constructor from a vector of colvars
2014-12-02 10:09:53 +08:00
colvar_grid ( std : : vector < colvar * > const & colvars ,
2015-02-20 02:20:52 +08:00
T const & t = T ( ) ,
size_t mult_i = 1 ,
bool margin = false )
: has_data ( false )
2012-05-24 00:20:04 +08:00
{
save_delimiters = false ;
2015-02-20 02:20:52 +08:00
this - > init_from_colvars ( colvars , t , mult_i , margin ) ;
}
int init_from_colvars ( std : : vector < colvar * > const & colvars ,
T const & t = T ( ) ,
size_t mult_i = 1 ,
bool margin = false )
{
if ( cvm : : debug ( ) ) {
cvm : : log ( " Reading grid configuration from collective variables. \n " ) ;
}
2012-05-24 00:20:04 +08:00
2015-02-20 02:20:52 +08:00
cv = colvars ;
nd = colvars . size ( ) ;
mult = mult_i ;
2012-05-24 00:20:04 +08:00
2015-02-20 02:20:52 +08:00
size_t i ;
if ( cvm : : debug ( ) ) {
2014-12-02 10:09:53 +08:00
cvm : : log ( " Allocating a grid for " + cvm : : to_str ( colvars . size ( ) ) +
2015-02-20 02:20:52 +08:00
" collective variables, multiplicity = " + cvm : : to_str ( mult_i ) + " . \n " ) ;
}
2012-05-24 00:20:04 +08:00
2015-02-20 02:20:52 +08:00
for ( i = 0 ; i < cv . size ( ) ; i + + ) {
2012-05-24 00:20:04 +08:00
2014-12-02 10:09:53 +08:00
if ( cv [ i ] - > value ( ) . type ( ) ! = colvarvalue : : type_scalar ) {
cvm : : error ( " Colvar grids can only be automatically "
2015-02-20 02:20:52 +08:00
" constructed for scalar variables. "
" ABF and histogram can not be used; metadynamics "
" can be used with useGrids disabled. \n " , INPUT_ERROR ) ;
return COLVARS_ERROR ;
2012-05-24 00:20:04 +08:00
}
if ( cv [ i ] - > width < = 0.0 ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Tried to initialize a grid on a "
2015-02-20 02:20:52 +08:00
" variable with negative or zero width. \n " , INPUT_ERROR ) ;
return COLVARS_ERROR ;
2012-05-24 00:20:04 +08:00
}
2014-12-02 10:09:53 +08:00
widths . push_back ( cv [ i ] - > width ) ;
hard_lower_boundaries . push_back ( cv [ i ] - > hard_lower_boundary ) ;
hard_upper_boundaries . push_back ( cv [ i ] - > hard_upper_boundary ) ;
periodic . push_back ( cv [ i ] - > periodic_boundaries ( ) ) ;
2012-05-24 00:20:04 +08:00
2012-06-30 01:52:31 +08:00
// By default, get reported colvar value (for extended Lagrangian colvars)
2014-12-02 10:09:53 +08:00
actual_value . push_back ( false ) ;
2012-06-30 01:52:31 +08:00
// except if a colvar is specified twice in a row
// then the first instance is the actual value
2016-08-04 03:32:03 +08:00
// For histograms of extended-system coordinates
2012-06-30 01:52:31 +08:00
if ( i > 0 & & cv [ i - 1 ] = = cv [ i ] ) {
actual_value [ i - 1 ] = true ;
}
2012-05-24 00:20:04 +08:00
if ( margin ) {
if ( periodic [ i ] ) {
// Shift the grid by half the bin width (values at edges instead of center of bins)
2014-12-02 10:09:53 +08:00
lower_boundaries . push_back ( cv [ i ] - > lower_boundary . real_value - 0.5 * widths [ i ] ) ;
upper_boundaries . push_back ( cv [ i ] - > upper_boundary . real_value - 0.5 * widths [ i ] ) ;
2012-05-24 00:20:04 +08:00
} else {
// Make this grid larger by one bin width
2014-12-02 10:09:53 +08:00
lower_boundaries . push_back ( cv [ i ] - > lower_boundary . real_value - 0.5 * widths [ i ] ) ;
upper_boundaries . push_back ( cv [ i ] - > upper_boundary . real_value + 0.5 * widths [ i ] ) ;
2012-05-24 00:20:04 +08:00
}
} else {
2014-12-02 10:09:53 +08:00
lower_boundaries . push_back ( cv [ i ] - > lower_boundary ) ;
upper_boundaries . push_back ( cv [ i ] - > upper_boundary ) ;
2012-05-24 00:20:04 +08:00
}
2015-02-20 02:20:52 +08:00
}
2012-05-24 00:20:04 +08:00
2015-02-20 02:20:52 +08:00
this - > init_from_boundaries ( ) ;
return this - > setup ( ) ;
}
2012-05-24 00:20:04 +08:00
2015-02-20 02:20:52 +08:00
int init_from_boundaries ( T const & t = T ( ) ,
size_t const & mult_i = 1 )
{
if ( cvm : : debug ( ) ) {
cvm : : log ( " Configuring grid dimensions from colvars boundaries. \n " ) ;
}
2013-06-28 06:48:27 +08:00
2015-02-20 02:20:52 +08:00
// these will have to be updated
nx . clear ( ) ;
nxc . clear ( ) ;
nt = 0 ;
for ( size_t i = 0 ; i < lower_boundaries . size ( ) ; i + + ) {
cvm : : real nbins = ( upper_boundaries [ i ] . real_value -
lower_boundaries [ i ] . real_value ) / widths [ i ] ;
int nbins_round = ( int ) ( nbins + 0.5 ) ;
if ( std : : fabs ( nbins - cvm : : real ( nbins_round ) ) > 1.0E-10 ) {
cvm : : log ( " Warning: grid interval( " +
cvm : : to_str ( lower_boundaries [ i ] , cvm : : cv_width , cvm : : cv_prec ) + " - " +
cvm : : to_str ( upper_boundaries [ i ] , cvm : : cv_width , cvm : : cv_prec ) +
" ) is not commensurate to its bin width( " +
cvm : : to_str ( widths [ i ] , cvm : : cv_width , cvm : : cv_prec ) + " ). \n " ) ;
upper_boundaries [ i ] . real_value = lower_boundaries [ i ] . real_value +
( nbins_round * widths [ i ] ) ;
2012-05-24 00:20:04 +08:00
}
2015-02-20 02:20:52 +08:00
if ( cvm : : debug ( ) )
cvm : : log ( " Number of points is " + cvm : : to_str ( ( int ) nbins_round ) +
" for the colvar no. " + cvm : : to_str ( i + 1 ) + " . \n " ) ;
nx . push_back ( nbins_round ) ;
2012-05-24 00:20:04 +08:00
}
2015-02-20 02:20:52 +08:00
return COLVARS_OK ;
2012-05-24 00:20:04 +08:00
}
/// Wrap an index vector around periodic boundary conditions
/// also checks validity of non-periodic indices
2014-12-02 10:09:53 +08:00
inline void wrap ( std : : vector < int > & ix ) const
2012-05-24 00:20:04 +08:00
{
for ( size_t i = 0 ; i < nd ; i + + ) {
if ( periodic [ i ] ) {
ix [ i ] = ( ix [ i ] + nx [ i ] ) % nx [ i ] ; //to ensure non-negative result
} else {
2015-02-20 02:20:52 +08:00
if ( ix [ i ] < 0 | | ix [ i ] > = nx [ i ] ) {
cvm : : error ( " Trying to wrap illegal index vector (non-PBC) for a grid point: "
+ cvm : : to_str ( ix ) , BUG_ERROR ) ;
2014-10-08 04:30:25 +08:00
return ;
2015-02-20 02:20:52 +08:00
}
2012-05-24 00:20:04 +08:00
}
}
}
2015-02-20 02:20:52 +08:00
2012-05-24 00:20:04 +08:00
/// \brief Report the bin corresponding to the current value of variable i
inline int current_bin_scalar ( int const i ) const
{
2014-12-02 10:09:53 +08:00
return value_to_bin_scalar ( actual_value [ i ] ? cv [ i ] - > actual_value ( ) : cv [ i ] - > value ( ) , i ) ;
2012-05-24 00:20:04 +08:00
}
2016-07-16 06:25:17 +08:00
/// \brief Report the bin corresponding to the current value of item iv in variable i
inline int current_bin_scalar ( int const i , int const iv ) const
{
return value_to_bin_scalar ( actual_value [ i ] ?
cv [ i ] - > actual_value ( ) . vector1d_value [ iv ] :
cv [ i ] - > value ( ) . vector1d_value [ iv ] , i ) ;
}
2012-05-24 00:20:04 +08:00
/// \brief Use the lower boundary and the width to report which bin
/// the provided value is in
2014-12-02 10:09:53 +08:00
inline int value_to_bin_scalar ( colvarvalue const & value , const int i ) const
2012-05-24 00:20:04 +08:00
{
2014-12-02 10:09:53 +08:00
return ( int ) std : : floor ( ( value . real_value - lower_boundaries [ i ] . real_value ) / widths [ i ] ) ;
2012-05-24 00:20:04 +08:00
}
/// \brief Same as the standard version, but uses another grid definition
2014-12-02 10:09:53 +08:00
inline int value_to_bin_scalar ( colvarvalue const & value ,
2012-05-24 00:20:04 +08:00
colvarvalue const & new_offset ,
cvm : : real const & new_width ) const
{
2014-12-02 10:09:53 +08:00
return ( int ) std : : floor ( ( value . real_value - new_offset . real_value ) / new_width ) ;
2012-05-24 00:20:04 +08:00
}
/// \brief Use the two boundaries and the width to report the
/// central value corresponding to a bin index
2014-12-02 10:09:53 +08:00
inline colvarvalue bin_to_value_scalar ( int const & i_bin , int const i ) const
2012-05-24 00:20:04 +08:00
{
return lower_boundaries [ i ] . real_value + widths [ i ] * ( 0.5 + i_bin ) ;
}
2013-06-28 06:48:27 +08:00
2012-05-24 00:20:04 +08:00
/// \brief Same as the standard version, but uses different parameters
2014-12-02 10:09:53 +08:00
inline colvarvalue bin_to_value_scalar ( int const & i_bin ,
2012-05-24 00:20:04 +08:00
colvarvalue const & new_offset ,
cvm : : real const & new_width ) const
{
return new_offset . real_value + new_width * ( 0.5 + i_bin ) ;
}
/// Set the value at the point with index ix
2014-12-02 10:09:53 +08:00
inline void set_value ( std : : vector < int > const & ix ,
2012-05-24 00:20:04 +08:00
T const & t ,
size_t const & imult = 0 )
2013-06-28 06:48:27 +08:00
{
2014-12-02 10:09:53 +08:00
data [ this - > address ( ix ) + imult ] = t ;
2012-05-24 00:20:04 +08:00
has_data = true ;
}
2014-10-29 03:53:17 +08:00
/// \brief Get the change from this to other_grid
/// and store the result in this.
/// this_grid := other_grid - this_grid
/// Grids must have the same dimensions.
2014-12-02 10:09:53 +08:00
void delta_grid ( colvar_grid < T > const & other_grid )
2014-10-29 03:53:17 +08:00
{
if ( other_grid . multiplicity ( ) ! = this - > multiplicity ( ) ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error: trying to subtract two grids with "
2014-10-29 03:53:17 +08:00
" different multiplicity. \n " ) ;
return ;
2015-02-20 02:20:52 +08:00
}
2014-10-29 03:53:17 +08:00
if ( other_grid . data . size ( ) ! = this - > data . size ( ) ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error: trying to subtract two grids with "
2014-10-29 03:53:17 +08:00
" different size. \n " ) ;
return ;
}
for ( size_t i = 0 ; i < data . size ( ) ; i + + ) {
data [ i ] = other_grid . data [ i ] - data [ i ] ;
}
has_data = true ;
}
/// \brief Copy data from another grid of the same type, AND
/// identical definition (boundaries, widths)
/// Added for shared ABF.
2014-12-02 10:09:53 +08:00
void copy_grid ( colvar_grid < T > const & other_grid )
2014-10-29 03:53:17 +08:00
{
if ( other_grid . multiplicity ( ) ! = this - > multiplicity ( ) ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error: trying to copy two grids with "
2014-10-29 03:53:17 +08:00
" different multiplicity. \n " ) ;
2015-02-20 02:20:52 +08:00
return ;
}
2014-10-29 03:53:17 +08:00
if ( other_grid . data . size ( ) ! = this - > data . size ( ) ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error: trying to copy two grids with "
2014-10-29 03:53:17 +08:00
" different size. \n " ) ;
return ;
}
for ( size_t i = 0 ; i < data . size ( ) ; i + + ) {
data [ i ] = other_grid . data [ i ] ;
}
has_data = true ;
}
/// \brief Extract the grid data as they are represented in memory.
/// Put the results in "out_data".
2014-12-02 10:09:53 +08:00
void raw_data_out ( T * out_data ) const
2014-10-29 03:53:17 +08:00
{
for ( size_t i = 0 ; i < data . size ( ) ; i + + ) out_data [ i ] = data [ i ] ;
}
/// \brief Input the data as they are represented in memory.
2014-12-02 10:09:53 +08:00
void raw_data_in ( const T * in_data )
2014-10-29 03:53:17 +08:00
{
for ( size_t i = 0 ; i < data . size ( ) ; i + + ) data [ i ] = in_data [ i ] ;
has_data = true ;
}
/// \brief Size of the data as they are represented in memory.
size_t raw_data_num ( ) const { return data . size ( ) ; }
2012-05-24 00:20:04 +08:00
/// \brief Get the binned value indexed by ix, or the first of them
/// if the multiplicity is larger than 1
2014-12-02 10:09:53 +08:00
inline T const & value ( std : : vector < int > const & ix ,
2012-05-24 00:20:04 +08:00
size_t const & imult = 0 ) const
{
2014-12-02 10:09:53 +08:00
return data [ this - > address ( ix ) + imult ] ;
2012-05-24 00:20:04 +08:00
}
/// \brief Add a constant to all elements (fast loop)
2014-12-02 10:09:53 +08:00
inline void add_constant ( T const & t )
2012-05-24 00:20:04 +08:00
{
2013-06-28 06:48:27 +08:00
for ( size_t i = 0 ; i < nt ; i + + )
2012-05-24 00:20:04 +08:00
data [ i ] + = t ;
has_data = true ;
}
/// \brief Multiply all elements by a scalar constant (fast loop)
2014-12-02 10:09:53 +08:00
inline void multiply_constant ( cvm : : real const & a )
2012-05-24 00:20:04 +08:00
{
2013-06-28 06:48:27 +08:00
for ( size_t i = 0 ; i < nt ; i + + )
2012-05-24 00:20:04 +08:00
data [ i ] * = a ;
}
2013-06-28 06:48:27 +08:00
2012-05-24 00:20:04 +08:00
/// \brief Get the bin indices corresponding to the provided values of
/// the colvars
2014-12-02 10:09:53 +08:00
inline std : : vector < int > const get_colvars_index ( std : : vector < colvarvalue > const & values ) const
2012-05-24 00:20:04 +08:00
{
std : : vector < int > index = new_index ( ) ;
for ( size_t i = 0 ; i < nd ; i + + ) {
2014-12-02 10:09:53 +08:00
index [ i ] = value_to_bin_scalar ( values [ i ] , i ) ;
2012-05-24 00:20:04 +08:00
}
return index ;
}
/// \brief Get the bin indices corresponding to the current values
/// of the colvars
inline std : : vector < int > const get_colvars_index ( ) const
{
std : : vector < int > index = new_index ( ) ;
for ( size_t i = 0 ; i < nd ; i + + ) {
2014-12-02 10:09:53 +08:00
index [ i ] = current_bin_scalar ( i ) ;
2012-05-24 00:20:04 +08:00
}
return index ;
}
/// \brief Get the minimal distance (in number of bins) from the
/// boundaries; a negative number is returned if the given point is
/// off-grid
2014-12-02 10:09:53 +08:00
inline cvm : : real bin_distance_from_boundaries ( std : : vector < colvarvalue > const & values ,
2012-06-30 01:52:31 +08:00
bool skip_hard_boundaries = false )
2012-05-24 00:20:04 +08:00
{
cvm : : real minimum = 1.0E+16 ;
for ( size_t i = 0 ; i < nd ; i + + ) {
if ( periodic [ i ] ) continue ;
2014-12-02 10:09:53 +08:00
cvm : : real dl = std : : sqrt ( cv [ i ] - > dist2 ( values [ i ] , lower_boundaries [ i ] ) ) / widths [ i ] ;
cvm : : real du = std : : sqrt ( cv [ i ] - > dist2 ( values [ i ] , upper_boundaries [ i ] ) ) / widths [ i ] ;
2012-05-24 00:20:04 +08:00
if ( values [ i ] . real_value < lower_boundaries [ i ] )
dl * = - 1.0 ;
if ( values [ i ] . real_value > upper_boundaries [ i ] )
du * = - 1.0 ;
2012-06-30 01:52:31 +08:00
if ( ( ( ! skip_hard_boundaries ) | | ( ! hard_lower_boundaries [ i ] ) ) & & ( dl < minimum ) )
2012-05-24 00:20:04 +08:00
minimum = dl ;
2012-06-30 01:52:31 +08:00
if ( ( ( ! skip_hard_boundaries ) | | ( ! hard_upper_boundaries [ i ] ) ) & & ( du < minimum ) )
2012-05-24 00:20:04 +08:00
minimum = du ;
}
return minimum ;
}
/// \brief Add data from another grid of the same type
2013-06-28 06:48:27 +08:00
///
2012-05-24 00:20:04 +08:00
/// Note: this function maps other_grid inside this one regardless
/// of whether it fits or not.
2014-12-02 10:09:53 +08:00
void map_grid ( colvar_grid < T > const & other_grid )
2012-05-24 00:20:04 +08:00
{
2014-10-08 04:30:25 +08:00
if ( other_grid . multiplicity ( ) ! = this - > multiplicity ( ) ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error: trying to merge two grids with values of "
2012-05-24 00:20:04 +08:00
" different multiplicity. \n " ) ;
2014-10-08 04:30:25 +08:00
return ;
}
2012-05-24 00:20:04 +08:00
std : : vector < colvarvalue > const & gb = this - > lower_boundaries ;
std : : vector < cvm : : real > const & gw = this - > widths ;
std : : vector < colvarvalue > const & ogb = other_grid . lower_boundaries ;
std : : vector < cvm : : real > const & ogw = other_grid . widths ;
std : : vector < int > ix = this - > new_index ( ) ;
std : : vector < int > oix = other_grid . new_index ( ) ;
if ( cvm : : debug ( ) )
2014-12-02 10:09:53 +08:00
cvm : : log ( " Remapping grid... \n " ) ;
for ( ; this - > index_ok ( ix ) ; this - > incr ( ix ) ) {
2012-05-24 00:20:04 +08:00
for ( size_t i = 0 ; i < nd ; i + + ) {
oix [ i ] =
2014-12-02 10:09:53 +08:00
value_to_bin_scalar ( bin_to_value_scalar ( ix [ i ] , gb [ i ] , gw [ i ] ) ,
2012-05-24 00:20:04 +08:00
ogb [ i ] ,
ogw [ i ] ) ;
}
2014-12-02 10:09:53 +08:00
if ( ! other_grid . index_ok ( oix ) ) {
2012-05-24 00:20:04 +08:00
continue ;
}
for ( size_t im = 0 ; im < mult ; im + + ) {
2014-12-02 10:09:53 +08:00
this - > set_value ( ix , other_grid . value ( oix , im ) , im ) ;
2012-05-24 00:20:04 +08:00
}
}
has_data = true ;
if ( cvm : : debug ( ) )
2014-12-02 10:09:53 +08:00
cvm : : log ( " Remapping done. \n " ) ;
2012-05-24 00:20:04 +08:00
}
/// \brief Add data from another grid of the same type, AND
/// identical definition (boundaries, widths)
2014-12-02 10:09:53 +08:00
void add_grid ( colvar_grid < T > const & other_grid ,
2012-05-24 00:20:04 +08:00
cvm : : real scale_factor = 1.0 )
{
2014-10-08 04:30:25 +08:00
if ( other_grid . multiplicity ( ) ! = this - > multiplicity ( ) ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error: trying to sum togetehr two grids with values of "
2012-05-24 00:20:04 +08:00
" different multiplicity. \n " ) ;
2014-10-08 04:30:25 +08:00
return ;
}
2013-06-28 06:48:27 +08:00
if ( scale_factor ! = 1.0 )
2012-05-24 00:20:04 +08:00
for ( size_t i = 0 ; i < data . size ( ) ; i + + ) {
data [ i ] + = scale_factor * other_grid . data [ i ] ;
}
2013-06-28 06:48:27 +08:00
else
2012-05-24 00:20:04 +08:00
// skip multiplication if possible
for ( size_t i = 0 ; i < data . size ( ) ; i + + ) {
data [ i ] + = other_grid . data [ i ] ;
}
has_data = true ;
}
/// \brief Return the value suitable for output purposes (so that it
/// may be rescaled or manipulated without changing it permanently)
2014-12-02 10:09:53 +08:00
virtual inline T value_output ( std : : vector < int > const & ix ,
2012-05-24 00:20:04 +08:00
size_t const & imult = 0 )
{
2014-12-02 10:09:53 +08:00
return value ( ix , imult ) ;
2012-05-24 00:20:04 +08:00
}
/// \brief 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)
2014-12-02 10:09:53 +08:00
virtual inline void value_input ( std : : vector < int > const & ix ,
2012-05-24 00:20:04 +08:00
T const & t ,
size_t const & imult = 0 ,
bool add = false )
{
if ( add )
2014-12-02 10:09:53 +08:00
data [ address ( ix ) + imult ] + = t ;
2012-05-24 00:20:04 +08:00
else
2014-12-02 10:09:53 +08:00
data [ address ( ix ) + imult ] = t ;
2012-05-24 00:20:04 +08:00
has_data = true ;
}
// /// Get the pointer to the binned value indexed by ix
// inline T const *value_p (std::vector<int> const &ix)
// {
// return &(data[address (ix)]);
// }
/// \brief Get the index corresponding to the "first" bin, to be
/// used as the initial value for an index in looping
inline std : : vector < int > const new_index ( ) const
{
return std : : vector < int > ( nd , 0 ) ;
}
/// \brief Check that the index is within range in each of the
/// dimensions
2014-12-02 10:09:53 +08:00
inline bool index_ok ( std : : vector < int > const & ix ) const
2012-05-24 00:20:04 +08:00
{
for ( size_t i = 0 ; i < nd ; i + + ) {
2014-12-02 10:09:53 +08:00
if ( ( ix [ i ] < 0 ) | | ( ix [ i ] > = int ( nx [ i ] ) ) )
2012-05-24 00:20:04 +08:00
return false ;
}
return true ;
}
/// \brief Increment the index, in a way that will make it loop over
/// the whole nd-dimensional array
2014-12-02 10:09:53 +08:00
inline void incr ( std : : vector < int > & ix ) const
2012-05-24 00:20:04 +08:00
{
for ( int i = ix . size ( ) - 1 ; i > = 0 ; i - - ) {
ix [ i ] + + ;
if ( ix [ i ] > = nx [ i ] ) {
if ( i > 0 ) {
ix [ i ] = 0 ;
continue ;
} else {
// this is the last iteration, a non-valid index is being
// set for the outer index, which will be caught by
// index_ok()
ix [ 0 ] = nx [ 0 ] ;
return ;
}
} else {
return ;
}
}
}
/// \brief Write the grid parameters (number of colvars, boundaries, width and number of points)
2014-12-02 10:09:53 +08:00
std : : ostream & write_params ( std : : ostream & os )
2012-05-24 00:20:04 +08:00
{
2014-10-08 04:30:25 +08:00
size_t i ;
2012-05-24 00:20:04 +08:00
os < < " grid_parameters { \n n_colvars " < < nd < < " \n " ;
2013-06-28 06:48:27 +08:00
2012-05-24 00:20:04 +08:00
os < < " lower_boundaries " ;
2014-10-08 04:30:25 +08:00
for ( i = 0 ; i < nd ; i + + )
2012-05-24 00:20:04 +08:00
os < < " " < < lower_boundaries [ i ] ;
os < < " \n " ;
os < < " upper_boundaries " ;
2014-10-08 04:30:25 +08:00
for ( i = 0 ; i < nd ; i + + )
2012-05-24 00:20:04 +08:00
os < < " " < < upper_boundaries [ i ] ;
os < < " \n " ;
os < < " widths " ;
2014-10-08 04:30:25 +08:00
for ( i = 0 ; i < nd ; i + + )
2012-05-24 00:20:04 +08:00
os < < " " < < widths [ i ] ;
os < < " \n " ;
os < < " sizes " ;
2014-10-08 04:30:25 +08:00
for ( i = 0 ; i < nd ; i + + )
2012-05-24 00:20:04 +08:00
os < < " " < < nx [ i ] ;
os < < " \n " ;
os < < " } \n " ;
return os ;
2013-06-28 06:48:27 +08:00
}
2012-05-24 00:20:04 +08:00
2015-02-20 02:20:52 +08:00
/// Read a grid definition from a config string
int parse_params ( std : : string const & conf )
2012-05-24 00:20:04 +08:00
{
2015-02-20 02:20:52 +08:00
if ( cvm : : debug ( ) ) cvm : : log ( " Reading grid configuration from string. \n " ) ;
2012-05-24 00:20:04 +08:00
std : : vector < int > old_nx = nx ;
std : : vector < colvarvalue > old_lb = lower_boundaries ;
{
size_t nd_in = 0 ;
2014-12-02 10:09:53 +08:00
colvarparse : : get_keyval ( conf , " n_colvars " , nd_in , nd , colvarparse : : parse_silent ) ;
2014-10-08 04:30:25 +08:00
if ( nd_in ! = nd ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error: trying to read data for a grid "
2014-10-08 04:30:25 +08:00
" that contains a different number of colvars ( " +
2014-12-02 10:09:53 +08:00
cvm : : to_str ( nd_in ) + " ) than the grid defined "
" in the configuration file( " + cvm : : to_str ( nd ) +
2014-10-08 04:30:25 +08:00
" ). \n " ) ;
2015-02-20 02:20:52 +08:00
return COLVARS_ERROR ;
2014-10-08 04:30:25 +08:00
}
2012-05-24 00:20:04 +08:00
}
2014-12-02 10:09:53 +08:00
colvarparse : : get_keyval ( conf , " lower_boundaries " ,
2012-05-24 00:20:04 +08:00
lower_boundaries , lower_boundaries , colvarparse : : parse_silent ) ;
2014-12-02 10:09:53 +08:00
colvarparse : : get_keyval ( conf , " upper_boundaries " ,
2012-05-24 00:20:04 +08:00
upper_boundaries , upper_boundaries , colvarparse : : parse_silent ) ;
2015-02-20 02:20:52 +08:00
// support also camel case
colvarparse : : get_keyval ( conf , " lowerBoundaries " ,
lower_boundaries , lower_boundaries , colvarparse : : parse_silent ) ;
colvarparse : : get_keyval ( conf , " upperBoundaries " ,
upper_boundaries , upper_boundaries , colvarparse : : parse_silent ) ;
2014-12-02 10:09:53 +08:00
colvarparse : : get_keyval ( conf , " widths " , widths , widths , colvarparse : : parse_silent ) ;
2012-05-24 00:20:04 +08:00
2014-12-02 10:09:53 +08:00
colvarparse : : get_keyval ( conf , " sizes " , nx , nx , colvarparse : : parse_silent ) ;
2012-05-24 00:20:04 +08:00
2015-02-20 02:20:52 +08:00
if ( nd < lower_boundaries . size ( ) ) nd = lower_boundaries . size ( ) ;
if ( ! actual_value . size ( ) ) actual_value . assign ( nd , false ) ;
if ( ! periodic . size ( ) ) periodic . assign ( nd , false ) ;
if ( ! widths . size ( ) ) widths . assign ( nd , 1.0 ) ;
2012-05-24 00:20:04 +08:00
bool new_params = false ;
2015-02-20 02:20:52 +08:00
if ( old_nx . size ( ) ) {
for ( size_t i = 0 ; i < nd ; i + + ) {
if ( ( old_nx [ i ] ! = nx [ i ] ) | |
( std : : sqrt ( cv [ i ] - > dist2 ( old_lb [ i ] ,
2012-05-24 00:20:04 +08:00
lower_boundaries [ i ] ) ) > 1.0E-10 ) ) {
2015-02-20 02:20:52 +08:00
new_params = true ;
}
2012-05-24 00:20:04 +08:00
}
2015-02-20 02:20:52 +08:00
} else {
new_params = true ;
2012-05-24 00:20:04 +08:00
}
// reallocate the array in case the grid params have just changed
if ( new_params ) {
2015-02-20 02:20:52 +08:00
init_from_boundaries ( ) ;
// data.resize(0); // no longer needed: setup calls clear()
return this - > setup ( nx , T ( ) , mult ) ;
2012-05-24 00:20:04 +08:00
}
2015-02-20 02:20:52 +08:00
return COLVARS_OK ;
2012-05-24 00:20:04 +08:00
}
/// \brief Check that the grid information inside (boundaries,
/// widths, ...) is consistent with the current setting of the
/// colvars
void check_consistency ( )
{
for ( size_t i = 0 ; i < nd ; i + + ) {
2014-12-02 10:09:53 +08:00
if ( ( std : : sqrt ( cv [ i ] - > dist2 ( cv [ i ] - > lower_boundary ,
2013-06-28 06:48:27 +08:00
lower_boundaries [ i ] ) ) > 1.0E-10 ) | |
2014-12-02 10:09:53 +08:00
( std : : sqrt ( cv [ i ] - > dist2 ( cv [ i ] - > upper_boundary ,
2013-06-28 06:48:27 +08:00
upper_boundaries [ i ] ) ) > 1.0E-10 ) | |
2014-12-02 10:09:53 +08:00
( std : : sqrt ( cv [ i ] - > dist2 ( cv [ i ] - > width ,
2012-05-24 00:20:04 +08:00
widths [ i ] ) ) > 1.0E-10 ) ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error: restart information for a grid is "
2014-10-08 04:30:25 +08:00
" inconsistent with that of its colvars. \n " ) ;
return ;
2012-05-24 00:20:04 +08:00
}
}
}
/// \brief Check that the grid information inside (boundaries,
2014-02-02 02:11:18 +08:00
/// widths, ...) is consistent with that of another grid
2014-12-02 10:09:53 +08:00
void check_consistency ( colvar_grid < T > const & other_grid )
2012-05-24 00:20:04 +08:00
{
for ( size_t i = 0 ; i < nd ; i + + ) {
// we skip dist2(), because periodicities and the like should
// matter: boundaries should be EXACTLY the same (otherwise,
// map_grid() should be used)
2014-12-02 10:09:53 +08:00
if ( ( std : : fabs ( other_grid . lower_boundaries [ i ] -
2013-06-28 06:48:27 +08:00
lower_boundaries [ i ] ) > 1.0E-10 ) | |
2014-12-02 10:09:53 +08:00
( std : : fabs ( other_grid . upper_boundaries [ i ] -
2013-06-28 06:48:27 +08:00
upper_boundaries [ i ] ) > 1.0E-10 ) | |
2014-12-02 10:09:53 +08:00
( std : : fabs ( other_grid . widths [ i ] -
2013-06-28 06:48:27 +08:00
widths [ i ] ) > 1.0E-10 ) | |
2012-05-24 00:20:04 +08:00
( data . size ( ) ! = other_grid . data . size ( ) ) ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error: inconsistency between "
2014-10-08 04:30:25 +08:00
" two grids that are supposed to be equal, "
" aside from the data stored. \n " ) ;
return ;
2012-05-24 00:20:04 +08:00
}
}
}
2013-06-28 06:48:27 +08:00
2012-05-24 00:20:04 +08:00
2015-09-09 07:38:11 +08:00
/// \brief Read grid entry in restart file
std : : istream & read_restart ( std : : istream & is )
{
size_t const start_pos = is . tellg ( ) ;
std : : string key , conf ;
if ( ( is > > key ) & & ( key = = std : : string ( " grid_parameters " ) ) ) {
is . seekg ( start_pos , std : : ios : : beg ) ;
is > > colvarparse : : read_block ( " grid_parameters " , conf ) ;
parse_params ( conf ) ;
} else {
cvm : : log ( " Grid parameters are missing in the restart file, using those from the configuration. \n " ) ;
is . seekg ( start_pos , std : : ios : : beg ) ;
}
read_raw ( is ) ;
return is ;
}
/// \brief Write grid entry in restart file
std : : ostream & write_restart ( std : : ostream & os )
{
write_params ( os ) ;
write_raw ( os ) ;
return os ;
}
2012-05-24 00:20:04 +08:00
/// \brief Write the grid data without labels, as they are
/// represented in memory
/// \param buf_size Number of values per line
2014-12-02 10:09:53 +08:00
std : : ostream & write_raw ( std : : ostream & os ,
2012-05-24 00:20:04 +08:00
size_t const buf_size = 3 )
{
std : : streamsize const w = os . width ( ) ;
std : : streamsize const p = os . precision ( ) ;
std : : vector < int > ix = new_index ( ) ;
size_t count = 0 ;
2014-12-02 10:09:53 +08:00
for ( ; index_ok ( ix ) ; incr ( ix ) ) {
2012-05-24 00:20:04 +08:00
for ( size_t imult = 0 ; imult < mult ; imult + + ) {
os < < " "
2014-12-02 10:09:53 +08:00
< < std : : setw ( w ) < < std : : setprecision ( p )
< < value_output ( ix , imult ) ;
2012-05-24 00:20:04 +08:00
if ( ( ( + + count ) % buf_size ) = = 0 )
os < < " \n " ;
}
}
// write a final newline only if buffer is not empty
if ( ( count % buf_size ) ! = 0 )
os < < " \n " ;
return os ;
}
2014-10-08 04:30:25 +08:00
/// \brief Read data written by colvar_grid::write_raw()
2014-12-02 10:09:53 +08:00
std : : istream & read_raw ( std : : istream & is )
2014-10-08 04:30:25 +08:00
{
size_t const start_pos = is . tellg ( ) ;
2013-06-28 06:48:27 +08:00
2014-12-02 10:09:53 +08:00
for ( std : : vector < int > ix = new_index ( ) ; index_ok ( ix ) ; incr ( ix ) ) {
2014-10-08 04:30:25 +08:00
for ( size_t imult = 0 ; imult < mult ; imult + + ) {
T new_value ;
if ( is > > new_value ) {
2014-12-02 10:09:53 +08:00
value_input ( ix , new_value , imult ) ;
2014-10-08 04:30:25 +08:00
} else {
is . clear ( ) ;
2014-12-02 10:09:53 +08:00
is . seekg ( start_pos , std : : ios : : beg ) ;
is . setstate ( std : : ios : : failbit ) ;
cvm : : 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 " ) ;
2014-10-08 04:30:25 +08:00
return is ;
}
2012-05-24 00:20:04 +08:00
}
}
2014-10-08 04:30:25 +08:00
has_data = true ;
return is ;
2012-05-24 00:20:04 +08:00
}
2014-10-08 04:30:25 +08:00
/// \brief Write the grid in a format which is both human readable
/// and suitable for visualization e.g. with gnuplot
2014-12-02 10:09:53 +08:00
void write_multicol ( std : : ostream & os )
2014-10-08 04:30:25 +08:00
{
std : : streamsize const w = os . width ( ) ;
std : : streamsize const p = os . precision ( ) ;
2012-05-24 00:20:04 +08:00
2014-10-08 04:30:25 +08:00
// Data in the header: nColvars, then for each
// xiMin, dXi, nPoints, periodic
2012-05-24 00:20:04 +08:00
2014-12-02 10:09:53 +08:00
os < < std : : setw ( 2 ) < < " # " < < nd < < " \n " ;
2012-05-24 00:20:04 +08:00
for ( size_t i = 0 ; i < nd ; i + + ) {
2014-10-08 04:30:25 +08:00
os < < " # "
2014-12-02 10:09:53 +08:00
< < std : : setw ( 10 ) < < lower_boundaries [ i ]
< < std : : setw ( 10 ) < < widths [ i ]
< < std : : setw ( 10 ) < < nx [ i ] < < " "
2014-10-08 04:30:25 +08:00
< < periodic [ i ] < < " \n " ;
2012-05-24 00:20:04 +08:00
}
2015-02-20 02:20:52 +08:00
2014-12-02 10:09:53 +08:00
for ( std : : vector < int > ix = new_index ( ) ; index_ok ( ix ) ; incr ( ix ) ) {
2012-05-24 00:20:04 +08:00
2014-10-08 04:30:25 +08:00
if ( ix . back ( ) = = 0 ) {
// if the last index is 0, add a new line to mark the new record
os < < " \n " ;
}
2012-05-24 00:20:04 +08:00
2014-10-08 04:30:25 +08:00
for ( size_t i = 0 ; i < nd ; i + + ) {
os < < " "
2014-12-02 10:09:53 +08:00
< < std : : setw ( w ) < < std : : setprecision ( p )
< < bin_to_value_scalar ( ix [ i ] , i ) ;
2014-10-08 04:30:25 +08:00
}
os < < " " ;
for ( size_t imult = 0 ; imult < mult ; imult + + ) {
os < < " "
2014-12-02 10:09:53 +08:00
< < std : : setw ( w ) < < std : : setprecision ( p )
< < value_output ( ix , imult ) ;
2014-10-08 04:30:25 +08:00
}
os < < " \n " ;
}
2012-05-24 00:20:04 +08:00
}
2014-10-08 04:30:25 +08:00
/// \brief Read a grid written by colvar_grid::write_multicol()
/// Adding data if add is true, replacing if false
2014-12-02 10:09:53 +08:00
std : : istream & read_multicol ( std : : istream & is , bool add = false )
2014-10-08 04:30:25 +08:00
{
// Data in the header: nColvars, then for each
// xiMin, dXi, nPoints, periodic
std : : string hash ;
cvm : : real lower , width , x ;
size_t n , periodic ;
bool remap ;
std : : vector < T > new_value ;
std : : vector < int > nx_read ;
std : : vector < int > bin ;
if ( cv . size ( ) ! = nd ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Cannot read grid file: missing reference to colvars. " ) ;
2014-10-08 04:30:25 +08:00
return is ;
}
2012-05-24 00:20:04 +08:00
if ( ! ( is > > hash ) | | ( hash ! = " # " ) ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error reading grid at position " +
cvm : : to_str ( is . tellg ( ) ) + " in stream(read \" " + hash + " \" ) \n " ) ;
2014-10-08 04:30:25 +08:00
return is ;
2012-05-24 00:20:04 +08:00
}
2014-10-08 04:30:25 +08:00
is > > n ;
if ( n ! = nd ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error reading grid: wrong number of collective variables. \n " ) ;
2014-10-08 04:30:25 +08:00
return is ;
}
2012-05-24 00:20:04 +08:00
2014-12-02 10:09:53 +08:00
nx_read . resize ( n ) ;
bin . resize ( n ) ;
new_value . resize ( mult ) ;
2012-05-24 00:20:04 +08:00
2014-10-08 04:30:25 +08:00
if ( this - > has_parent_data & & add ) {
2014-12-02 10:09:53 +08:00
new_data . resize ( data . size ( ) ) ;
2012-05-24 00:20:04 +08:00
}
2014-10-08 04:30:25 +08:00
remap = false ;
for ( size_t i = 0 ; i < nd ; i + + ) {
if ( ! ( is > > hash ) | | ( hash ! = " # " ) ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error reading grid at position " +
cvm : : to_str ( is . tellg ( ) ) + " in stream(read \" " + hash + " \" ) \n " ) ;
2014-10-08 04:30:25 +08:00
return is ;
2012-05-24 00:20:04 +08:00
}
2014-10-08 04:30:25 +08:00
is > > lower > > width > > nx_read [ i ] > > periodic ;
2014-12-02 10:09:53 +08:00
if ( ( std : : fabs ( lower - lower_boundaries [ i ] . real_value ) > 1.0e-10 ) | |
( std : : fabs ( width - widths [ i ] ) > 1.0e-10 ) | |
2014-10-08 04:30:25 +08:00
( nx_read [ i ] ! = nx [ i ] ) ) {
2014-12-02 10:09:53 +08:00
cvm : : log ( " Warning: reading from different grid definition (colvar "
+ cvm : : to_str ( i + 1 ) + " ); remapping data on new grid. \n " ) ;
2014-10-08 04:30:25 +08:00
remap = true ;
2012-05-24 00:20:04 +08:00
}
2014-10-08 04:30:25 +08:00
}
if ( remap ) {
// re-grid data
while ( is . good ( ) ) {
bool end_of_file = false ;
for ( size_t i = 0 ; i < nd ; i + + ) {
if ( ! ( is > > x ) ) end_of_file = true ;
2014-12-02 10:09:53 +08:00
bin [ i ] = value_to_bin_scalar ( x , i ) ;
2014-10-08 04:30:25 +08:00
}
if ( end_of_file ) break ;
2012-05-24 00:20:04 +08:00
for ( size_t imult = 0 ; imult < mult ; imult + + ) {
2014-10-08 04:30:25 +08:00
is > > new_value [ imult ] ;
}
if ( index_ok ( bin ) ) {
for ( size_t imult = 0 ; imult < mult ; imult + + ) {
2014-12-02 10:09:53 +08:00
value_input ( bin , new_value [ imult ] , imult , add ) ;
2014-10-08 04:30:25 +08:00
}
2012-05-24 00:20:04 +08:00
}
}
2014-10-08 04:30:25 +08:00
} else {
// do not re-grid the data but assume the same grid is used
2014-12-02 10:09:53 +08:00
for ( std : : vector < int > ix = new_index ( ) ; index_ok ( ix ) ; incr ( ix ) ) {
2014-10-08 04:30:25 +08:00
for ( size_t i = 0 ; i < nd ; i + + ) {
is > > x ;
}
for ( size_t imult = 0 ; imult < mult ; imult + + ) {
is > > new_value [ imult ] ;
2014-12-02 10:09:53 +08:00
value_input ( ix , new_value [ imult ] , imult , add ) ;
2014-10-08 04:30:25 +08:00
}
2012-05-24 00:20:04 +08:00
}
}
2014-10-08 04:30:25 +08:00
has_data = true ;
return is ;
2012-05-24 00:20:04 +08:00
}
2015-07-22 22:36:59 +08:00
/// \brief Write the grid data without labels, as they are
/// represented in memory
/// \param buf_size Number of values per line
std : : ostream & write_opendx ( std : : ostream & os )
{
// write the header
os < < " object 1 class gridpositions counts " ;
int icv ;
for ( icv = 0 ; icv < number_of_colvars ( ) ; icv + + ) {
os < < " " < < number_of_points ( icv ) ;
}
os < < " \n " ;
os < < " origin " ;
for ( icv = 0 ; icv < number_of_colvars ( ) ; icv + + ) {
os < < " " < < ( lower_boundaries [ icv ] . real_value + 0.5 * widths [ icv ] ) ;
}
os < < " \n " ;
for ( icv = 0 ; icv < number_of_colvars ( ) ; icv + + ) {
os < < " delta " ;
for ( size_t icv2 = 0 ; icv2 < number_of_colvars ( ) ; icv2 + + ) {
if ( icv = = icv2 ) os < < " " < < widths [ icv ] ;
else os < < " " < < 0.0 ;
}
os < < " \n " ;
}
os < < " object 2 class gridconnections counts " ;
for ( icv = 0 ; icv < number_of_colvars ( ) ; icv + + ) {
os < < " " < < number_of_points ( icv ) ;
}
os < < " \n " ;
os < < " object 3 class array type double rank 0 items "
< < number_of_points ( ) < < " data follows \n " ;
write_raw ( os ) ;
os < < " object \" collective variables scalar field \" class field \n " ;
return os ;
}
2012-05-24 00:20:04 +08:00
} ;
/// \brief Colvar_grid derived class to hold counters in discrete
/// n-dim colvar space
class colvar_grid_count : public colvar_grid < size_t >
{
public :
/// Default constructor
colvar_grid_count ( ) ;
/// Destructor
virtual inline ~ colvar_grid_count ( )
{ }
/// Constructor
2014-12-02 10:09:53 +08:00
colvar_grid_count ( std : : vector < int > const & nx_i ,
2012-05-24 00:20:04 +08:00
size_t const & def_count = 0 ) ;
/// Constructor from a vector of colvars
2014-12-02 10:09:53 +08:00
colvar_grid_count ( std : : vector < colvar * > & colvars ,
2012-05-24 00:20:04 +08:00
size_t const & def_count = 0 ) ;
/// Increment the counter at given position
2014-12-02 10:09:53 +08:00
inline void incr_count ( std : : vector < int > const & ix )
2012-05-24 00:20:04 +08:00
{
2014-12-02 10:09:53 +08:00
+ + ( data [ this - > address ( ix ) ] ) ;
2012-05-24 00:20:04 +08:00
}
2013-06-28 06:48:27 +08:00
/// \brief Get the binned count indexed by ix from the newly read data
2014-12-02 10:09:53 +08:00
inline size_t const & new_count ( std : : vector < int > const & ix ,
2012-05-24 00:20:04 +08:00
size_t const & imult = 0 )
{
2014-12-02 10:09:53 +08:00
return new_data [ address ( ix ) + imult ] ;
2012-05-24 00:20:04 +08:00
}
/// \brief Get the value from a formatted output and transform it
/// into the internal representation (it may have been rescaled or
/// manipulated)
2014-12-02 10:09:53 +08:00
virtual inline void value_input ( std : : vector < int > const & ix ,
2012-05-24 00:20:04 +08:00
size_t const & t ,
size_t const & imult = 0 ,
bool add = false )
{
if ( add ) {
2014-12-02 10:09:53 +08:00
data [ address ( ix ) ] + = t ;
2012-05-24 00:20:04 +08:00
if ( this - > has_parent_data ) {
// save newly read data for inputting parent grid
2014-12-02 10:09:53 +08:00
new_data [ address ( ix ) ] = t ;
2012-05-24 00:20:04 +08:00
}
} else {
2014-12-02 10:09:53 +08:00
data [ address ( ix ) ] = t ;
2012-05-24 00:20:04 +08:00
}
has_data = true ;
}
2016-09-07 09:26:50 +08:00
/// \brief Return the log-gradient from finite differences
/// on the *same* grid for dimension n
inline const cvm : : real log_gradient_finite_diff ( const std : : vector < int > & ix0 ,
int n = 0 )
{
cvm : : real A0 , A1 ;
std : : vector < int > ix ;
// factor for mesh width, 2.0 for central finite difference
// but only 1.0 on edges for non-PBC coordinates
cvm : : real factor ;
if ( periodic [ n ] ) {
factor = 2. ;
ix = ix0 ;
ix [ n ] - - ; wrap ( ix ) ;
A0 = data [ address ( ix ) ] ;
ix = ix0 ;
ix [ n ] + + ; wrap ( ix ) ;
A1 = data [ address ( ix ) ] ;
} else {
factor = 0. ;
ix = ix0 ;
if ( ix [ n ] > 0 ) { // not left edge
ix [ n ] - - ;
factor + = 1. ;
}
A0 = data [ address ( ix ) ] ;
ix = ix0 ;
if ( ix [ n ] + 1 < nx [ n ] ) { // not right edge
ix [ n ] + + ;
factor + = 1. ;
}
A1 = data [ address ( ix ) ] ;
}
if ( A0 = = 0 | | A1 = = 0 ) {
// can't handle empty bins
return 0. ;
} else {
return ( std : : log ( ( cvm : : real ) A1 ) - std : : log ( ( cvm : : real ) A0 ) )
/ ( widths [ n ] * factor ) ;
}
}
2012-05-24 00:20:04 +08:00
} ;
/// Class for accumulating a scalar function on a grid
class colvar_grid_scalar : public colvar_grid < cvm : : real >
{
public :
/// \brief Provide the associated sample count by which each binned value
/// should be divided
colvar_grid_count * samples ;
/// Default constructor
colvar_grid_scalar ( ) ;
/// Copy constructor (needed because of the grad pointer)
2014-12-02 10:09:53 +08:00
colvar_grid_scalar ( colvar_grid_scalar const & g ) ;
2012-05-24 00:20:04 +08:00
/// Destructor
~ colvar_grid_scalar ( ) ;
/// Constructor from specific sizes arrays
2014-12-02 10:09:53 +08:00
colvar_grid_scalar ( std : : vector < int > const & nx_i ) ;
2012-05-24 00:20:04 +08:00
/// Constructor from a vector of colvars
2014-12-02 10:09:53 +08:00
colvar_grid_scalar ( std : : vector < colvar * > & colvars ,
2015-07-22 22:36:59 +08:00
bool margin = 0 ) ;
2012-05-24 00:20:04 +08:00
/// Accumulate the value
2014-12-02 10:09:53 +08:00
inline void acc_value ( std : : vector < int > const & ix ,
2015-07-22 22:36:59 +08:00
cvm : : real const & new_value ,
size_t const & imult = 0 )
2012-05-24 00:20:04 +08:00
{
// only legal value of imult here is 0
2014-12-02 10:09:53 +08:00
data [ address ( ix ) ] + = new_value ;
2012-05-24 00:20:04 +08:00
if ( samples )
2014-12-02 10:09:53 +08:00
samples - > incr_count ( ix ) ;
2012-05-24 00:20:04 +08:00
has_data = true ;
}
/// Return the gradient of the scalar field from finite differences
2014-12-02 10:09:53 +08:00
inline const cvm : : real * gradient_finite_diff ( const std : : vector < int > & ix0 )
2012-05-24 00:20:04 +08:00
{
cvm : : real A0 , A1 ;
std : : vector < int > ix ;
2014-10-08 04:30:25 +08:00
if ( nd ! = 2 ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Finite differences available in dimension 2 only. " ) ;
2014-10-08 04:30:25 +08:00
return grad ;
}
for ( unsigned int n = 0 ; n < nd ; n + + ) {
2012-05-24 00:20:04 +08:00
ix = ix0 ;
2014-12-02 10:09:53 +08:00
A0 = data [ address ( ix ) ] ;
ix [ n ] + + ; wrap ( ix ) ;
A1 = data [ address ( ix ) ] ;
ix [ 1 - n ] + + ; wrap ( ix ) ;
A1 + = data [ address ( ix ) ] ;
ix [ n ] - - ; wrap ( ix ) ;
A0 + = data [ address ( ix ) ] ;
2012-05-24 00:20:04 +08:00
grad [ n ] = 0.5 * ( A1 - A0 ) / widths [ n ] ;
}
2013-06-28 06:48:27 +08:00
return grad ;
2012-05-24 00:20:04 +08:00
}
/// \brief Return the value of the function at ix divided by its
/// number of samples (if the count grid is defined)
2014-12-02 10:09:53 +08:00
virtual cvm : : real value_output ( std : : vector < int > const & ix ,
2015-07-22 22:36:59 +08:00
size_t const & imult = 0 )
2012-05-24 00:20:04 +08:00
{
2014-10-08 04:30:25 +08:00
if ( imult > 0 ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error: trying to access a component "
2015-07-22 22:36:59 +08:00
" larger than 1 in a scalar data grid. \n " ) ;
2014-10-08 04:30:25 +08:00
return 0. ;
}
2015-07-22 22:36:59 +08:00
if ( samples ) {
2014-12-02 10:09:53 +08:00
return ( samples - > value ( ix ) > 0 ) ?
( data [ address ( ix ) ] / cvm : : real ( samples - > value ( ix ) ) ) :
2012-05-24 00:20:04 +08:00
0.0 ;
2015-07-22 22:36:59 +08:00
} else {
2014-12-02 10:09:53 +08:00
return data [ address ( ix ) ] ;
2015-07-22 22:36:59 +08:00
}
2012-05-24 00:20:04 +08:00
}
/// \brief Get the value from a formatted output and transform it
/// into the internal representation (it may have been rescaled or
/// manipulated)
2014-12-02 10:09:53 +08:00
virtual void value_input ( std : : vector < int > const & ix ,
2015-07-22 22:36:59 +08:00
cvm : : real const & new_value ,
size_t const & imult = 0 ,
bool add = false )
2012-05-24 00:20:04 +08:00
{
2014-10-08 04:30:25 +08:00
if ( imult > 0 ) {
2014-12-02 10:09:53 +08:00
cvm : : error ( " Error: trying to access a component "
2015-07-22 22:36:59 +08:00
" larger than 1 in a scalar data grid. \n " ) ;
2014-10-08 04:30:25 +08:00
return ;
}
2012-05-24 00:20:04 +08:00
if ( add ) {
if ( samples )
2014-12-02 10:09:53 +08:00
data [ address ( ix ) ] + = new_value * samples - > new_count ( ix ) ;
2012-05-24 00:20:04 +08:00
else
2014-12-02 10:09:53 +08:00
data [ address ( ix ) ] + = new_value ;
2012-05-24 00:20:04 +08:00
} else {
if ( samples )
2014-12-02 10:09:53 +08:00
data [ address ( ix ) ] = new_value * samples - > value ( ix ) ;
2012-05-24 00:20:04 +08:00
else
2014-12-02 10:09:53 +08:00
data [ address ( ix ) ] = new_value ;
2012-05-24 00:20:04 +08:00
}
has_data = true ;
}
/// \brief Return the highest value
2015-07-22 22:36:59 +08:00
cvm : : real maximum_value ( ) const ;
2012-05-24 00:20:04 +08:00
/// \brief Return the lowest value
2015-07-22 22:36:59 +08:00
cvm : : real minimum_value ( ) const ;
/// \brief Calculates the integral of the map (uses widths if they are defined)
cvm : : real integral ( ) const ;
/// \brief Assuming that the map is a normalized probability density,
/// calculates the entropy (uses widths if they are defined)
cvm : : real entropy ( ) const ;
2012-05-24 00:20:04 +08:00
private :
// gradient
cvm : : real * grad ;
} ;
/// Class for accumulating the gradient of a scalar function on a grid
class colvar_grid_gradient : public colvar_grid < cvm : : real >
{
public :
/// \brief Provide the sample count by which each binned value
/// should be divided
colvar_grid_count * samples ;
/// Default constructor
colvar_grid_gradient ( ) ;
/// Destructor
virtual inline ~ colvar_grid_gradient ( )
{ }
/// Constructor from specific sizes arrays
2014-12-02 10:09:53 +08:00
colvar_grid_gradient ( std : : vector < int > const & nx_i ) ;
2012-05-24 00:20:04 +08:00
/// Constructor from a vector of colvars
2014-12-02 10:09:53 +08:00
colvar_grid_gradient ( std : : vector < colvar * > & colvars ) ;
2012-05-24 00:20:04 +08:00
/// \brief Accumulate the gradient
2014-12-02 10:09:53 +08:00
inline void acc_grad ( std : : vector < int > const & ix , cvm : : real const * grads ) {
2012-05-24 00:20:04 +08:00
for ( size_t imult = 0 ; imult < mult ; imult + + ) {
2014-12-02 10:09:53 +08:00
data [ address ( ix ) + imult ] + = grads [ imult ] ;
2012-05-24 00:20:04 +08:00
}
if ( samples )
2014-12-02 10:09:53 +08:00
samples - > incr_count ( ix ) ;
2012-05-24 00:20:04 +08:00
}
/// \brief Accumulate the gradient based on the force (i.e. sums the
/// opposite of the force)
2014-12-02 10:09:53 +08:00
inline void acc_force ( std : : vector < int > const & ix , cvm : : real const * forces ) {
2012-05-24 00:20:04 +08:00
for ( size_t imult = 0 ; imult < mult ; imult + + ) {
2014-12-02 10:09:53 +08:00
data [ address ( ix ) + imult ] - = forces [ imult ] ;
2012-05-24 00:20:04 +08:00
}
if ( samples )
2014-12-02 10:09:53 +08:00
samples - > incr_count ( ix ) ;
2012-05-24 00:20:04 +08:00
}
/// \brief Return the value of the function at ix divided by its
/// number of samples (if the count grid is defined)
2014-12-02 10:09:53 +08:00
virtual inline cvm : : real value_output ( std : : vector < int > const & ix ,
2012-05-24 00:20:04 +08:00
size_t const & imult = 0 )
{
if ( samples )
2014-12-02 10:09:53 +08:00
return ( samples - > value ( ix ) > 0 ) ?
( data [ address ( ix ) + imult ] / cvm : : real ( samples - > value ( ix ) ) ) :
2012-05-24 00:20:04 +08:00
0.0 ;
else
2014-12-02 10:09:53 +08:00
return data [ address ( ix ) + imult ] ;
2012-05-24 00:20:04 +08:00
}
/// \brief Get the value from a formatted output and transform it
/// into the internal representation (it may have been rescaled or
/// manipulated)
2014-12-02 10:09:53 +08:00
virtual inline void value_input ( std : : vector < int > const & ix ,
2012-05-24 00:20:04 +08:00
cvm : : real const & new_value ,
size_t const & imult = 0 ,
bool add = false )
{
if ( add ) {
if ( samples )
2014-12-02 10:09:53 +08:00
data [ address ( ix ) + imult ] + = new_value * samples - > new_count ( ix ) ;
2012-05-24 00:20:04 +08:00
else
2014-12-02 10:09:53 +08:00
data [ address ( ix ) + imult ] + = new_value ;
2012-05-24 00:20:04 +08:00
} else {
if ( samples )
2014-12-02 10:09:53 +08:00
data [ address ( ix ) + imult ] = new_value * samples - > value ( ix ) ;
2012-05-24 00:20:04 +08:00
else
2014-12-02 10:09:53 +08:00
data [ address ( ix ) + imult ] = new_value ;
2012-05-24 00:20:04 +08:00
}
has_data = true ;
}
/// Compute and return average value for a 1D gradient grid
inline cvm : : real average ( )
{
size_t n = 0 ;
if ( nd ! = 1 | | nx [ 0 ] = = 0 ) {
return 0.0 ;
}
cvm : : real sum = 0.0 ;
std : : vector < int > ix = new_index ( ) ;
if ( samples ) {
2014-12-02 10:09:53 +08:00
for ( ; index_ok ( ix ) ; incr ( ix ) ) {
if ( ( n = samples - > value ( ix ) ) )
sum + = value ( ix ) / n ;
2012-05-24 00:20:04 +08:00
}
} else {
2014-12-02 10:09:53 +08:00
for ( ; index_ok ( ix ) ; incr ( ix ) ) {
sum + = value ( ix ) ;
2012-05-24 00:20:04 +08:00
}
}
2014-12-02 10:09:53 +08:00
return ( sum / cvm : : real ( nx [ 0 ] ) ) ;
2012-05-24 00:20:04 +08:00
}
/// \brief If the grid is 1-dimensional, integrate it and write the
/// integral to a file
2014-12-02 10:09:53 +08:00
void write_1D_integral ( std : : ostream & os ) ;
2012-05-24 00:20:04 +08:00
} ;
# endif