Collective Variables Module - Developer Documentation
Loading...
Searching...
No Matches
colvargrid.h
1// -*- c++ -*-
2
3// This file is part of the Collective Variables module (Colvars).
4// The original version of Colvars and its updates are located at:
5// https://github.com/Colvars/colvars
6// Please update all Colvars source files before making any changes.
7// If you wish to distribute your changes, please submit them to the
8// Colvars repository at GitHub.
9
10#ifndef COLVARGRID_H
11#define COLVARGRID_H
12
13#include <iosfwd>
14#include <memory>
15
16#include "colvar.h"
17#include "colvarmodule.h"
18#include "colvarvalue.h"
19#include "colvarparse.h"
20
21
25
26public:
28 size_t nd = 0;
29
31 std::vector<int> nx;
32
34 std::vector<int> nxc;
35
37 std::vector<colvarvalue> lower_boundaries;
38
40 std::vector<colvarvalue> upper_boundaries;
41
43 std::vector<cvm::real> widths;
44};
45
46
52template <class T> class colvar_grid : public colvar_grid_params, public colvarparse {
53
54 //protected:
55public: // TODO create accessors for these after all instantiations work
56
59 size_t mult = 1;
60
62 size_t nt;
63
65 std::vector<T> data;
66
68 std::vector<size_t> new_data;
69
71 std::vector<colvar *> cv;
72
74 std::vector<bool> use_actual_value;
75
77 inline size_t address(std::vector<int> const &ix) const
78 {
79 size_t addr = 0;
80 for (size_t i = 0; i < nd; i++) {
81 addr += ix[i]*static_cast<size_t>(nxc[i]);
82 if (cvm::debug()) {
83 if (ix[i] >= nx[i]) {
84 cvmodule->error("Error: exceeding bounds in colvar_grid.\n", COLVARS_BUG_ERROR);
85 return 0;
86 }
87 }
88 }
89 return addr;
90 }
91
92public:
94 std::vector<bool> periodic;
95
97 std::vector<bool> hard_lower_boundaries;
98
100 std::vector<bool> hard_upper_boundaries;
101
104
107
109 inline size_t num_variables() const
110 {
111 return nd;
112 }
113
115 inline std::vector<int> const &number_of_points_vec() const
116 {
117 return nx;
118 }
119
122 inline size_t number_of_points(int const icv = -1) const
123 {
124 if (icv < 0) {
125 return nt;
126 } else {
127 return nx[icv];
128 }
129 }
130
132 inline std::vector<int> const & sizes() const
133 {
134 return nx;
135 }
136
138 inline void set_sizes(std::vector<int> const &new_sizes)
139 {
140 nx = new_sizes;
141 }
142
144 inline size_t multiplicity() const
145 {
146 return mult;
147 }
148
150 inline void request_actual_value(bool b = true)
151 {
152 size_t i;
153 for (i = 0; i < use_actual_value.size(); i++) {
154 use_actual_value[i] = b;
155 }
156 }
157
159 int setup(std::vector<int> const &nx_i,
160 T const &t = T(),
161 size_t const &mult_i = 1)
162 {
163 if (cvm::debug()) {
164 cvmodule->log("Allocating grid: multiplicity = "+cvm::to_str(mult_i)+
165 ", dimensions = "+cvm::to_str(nx_i)+".\n");
166 }
167
168 mult = mult_i;
169
170 data.clear();
171
172 nx = nx_i;
173 nd = nx.size();
174
175 nxc.resize(nd);
176
177 // setup dimensions
178 nt = mult;
179 for (int i = nd-1; i >= 0; i--) {
180 if (nx[i] <= 0) {
181 cvm::error_static("Error: providing an invalid number of grid points, "+
182 cvm::to_str(nx[i])+".\n", COLVARS_BUG_ERROR);
183 return COLVARS_ERROR;
184 }
185 nxc[i] = nt;
186 nt *= nx[i];
187 }
188
189 if (cvm::debug()) {
190 cvmodule->log("Total number of grid elements = "+cvm::to_str(nt)+".\n");
191 }
192
193 data.reserve(nt);
194 data.assign(nt, t);
195
196 return COLVARS_OK;
197 }
198
200 int setup()
201 {
202 return setup(this->nx, T(), this->mult);
203 }
204
206 void reset(T const &t = T())
207 {
208 data.assign(nt, t);
209 }
210
211
213 // This constructor depends on a static cvm pointer and is deprecated
214 colvar_grid() : colvarparse(cvm::main()), has_data(false)
215 {
216 nd = nt = 0;
217 mult = 1;
218 has_parent_data = false;
219 this->setup();
220 }
221
223 virtual ~colvar_grid()
224 {}
225
229 // This constructor depends on a static cvm pointer and is deprecated
231 colvarparse(cvm::main()),
232 mult(g.mult),
233 data(),
234 cv(g.cv),
239 has_parent_data(false),
240 has_data(false)
241 {}
242
247 colvar_grid(std::vector<int> const &nx_i,
248 T const &t = T(),
249 size_t mult_i = 1)
250 : colvarparse(cvm::main()), has_parent_data(false), has_data(false)
251 {
252 this->setup(nx_i, t, mult_i);
253 }
254
258 colvar_grid(std::vector<colvar *> const &colvars,
259 T const &t = T(),
260 size_t mult_i = 1,
261 bool add_extra_bin = false,
262 std::shared_ptr<const colvar_grid_params> params = nullptr,
263 std::string config = std::string())
264 : colvarparse(cvm::main()), has_parent_data(false), has_data(false)
265 {
266 (void) t;
267 this->init_from_colvars(colvars, mult_i, add_extra_bin, params, config);
268 }
269
273 colvar_grid(std::string const &filename, size_t mult_i = 1);
274
275 int init_from_colvars(std::vector<colvar *> const &colvars,
276 size_t mult_i = 1,
277 bool add_extra_bin = false,
278 std::shared_ptr<const colvar_grid_params> params = nullptr,
279 std::string config = std::string())
280 {
281 if (cvm::debug()) {
282 cvmodule->log("Reading grid configuration from collective variables.\n");
283 }
284
285 cv = colvars;
286 nd = colvars.size();
287 mult = mult_i;
288
289 size_t i;
290
291 if (cvm::debug()) {
292 cvmodule->log("Allocating a grid for "+cvm::to_str(colvars.size())+
293 " collective variables, multiplicity = "+cvm::to_str(mult_i)+".\n");
294 }
295
296 for (i = 0; i < nd; i++) {
297 if (cv[i]->value().type() != colvarvalue::type_scalar) {
298 cvm::error_static("Colvar grids can only be automatically "
299 "constructed for scalar variables. "
300 "ABF and histogram can not be used; metadynamics "
301 "can be used with useGrids disabled.\n", COLVARS_INPUT_ERROR);
302 return COLVARS_ERROR;
303 }
304
305 if (cv[i]->width <= 0.0) {
306 cvm::error_static("Tried to initialize a grid on a "
307 "variable with negative or zero width.\n", COLVARS_INPUT_ERROR);
308 return COLVARS_ERROR;
309 }
310
311 widths.push_back(cv[i]->width);
314
315 // By default, get reported colvar value (for extended Lagrangian colvars)
316 use_actual_value.push_back(false);
317
318 // except if a colvar is specified twice in a row
319 // then the first instance is the actual value
320 // For histograms of extended-system coordinates
321 if (i > 0 && cv[i-1] == cv[i]) {
322 use_actual_value[i-1] = true;
323 }
324
325 // This needs to work if the boundaries are undefined in the colvars
326 lower_boundaries.push_back(cv[i]->lower_boundary);
327 upper_boundaries.push_back(cv[i]->upper_boundary);
328 }
329
330 // Replace widths and boundaries with optional custom configuration
331 if (!config.empty()) {
332 this->parse_params(config);
333 this->check_keywords(config, "grid");
334
335 if (params) {
336 cvm::error_static("Error: init_from_colvars was passed both a grid config and a template grid.", COLVARS_BUG_ERROR);
337 return COLVARS_BUG_ERROR;
338 }
339 } else if (params) {
340 // Match grid sizes with template
341
342 if (params->nd != nd) {
343 cvm::error_static("Trying to initialize grid from template with wrong dimension (" +
344 cvm::to_str(params->nd) + " instead of " +
345 cvm::to_str(this->nd) + ").");
346 return COLVARS_ERROR;
347 }
348
349 widths =params->widths;
350 lower_boundaries =params->lower_boundaries;
351 upper_boundaries =params->upper_boundaries;
352 }
353
354 // Only now can we determine periodicity
355 for (i = 0; i < nd; i++) {
356 periodic.push_back(cv[i]->periodic_boundaries(lower_boundaries[i].real_value,
357 upper_boundaries[i].real_value));
358
359 if (add_extra_bin) {
360 // Shift the grid by half the bin width (values at edges instead of center of bins)
361 lower_boundaries[i] -= 0.5 * widths[i];
362
363 if (periodic[i]) {
364 // Just shift
365 upper_boundaries[i] -= 0.5 * widths[i];
366 } else {
367 // Widen grid by one bin width
368 upper_boundaries[i] += 0.5 * widths[i];
369 }
370 }
371 }
372
373 // Reset grid sizes based on widths and boundaries
374 this->init_from_boundaries();
375 return this->setup();
376 }
377
378 int init_from_boundaries()
379 {
380 if (cvm::debug()) {
381 cvmodule->log("Configuring grid dimensions from colvars boundaries.\n");
382 }
383
384 // these will have to be updated
385 nx.clear();
386 nxc.clear();
387 nt = 0;
388
389 for (size_t i = 0; i < lower_boundaries.size(); i++) {
390 // Re-compute periodicity using current grid boundaries
391 periodic[i] = cv[i]->periodic_boundaries(lower_boundaries[i].real_value,
392 upper_boundaries[i].real_value);
393
394 cvm::real nbins = ( upper_boundaries[i].real_value -
395 lower_boundaries[i].real_value ) / widths[i];
396 int nbins_round = (int)(nbins+0.5);
397
398 if (cvm::fabs(nbins - cvm::real(nbins_round)) > 1.0E-10) {
399 cvmodule->log("Warning: grid interval("+
400 cvm::to_str(lower_boundaries[i], cvmodule->cv_width, cvmodule->cv_prec)+" - "+
401 cvm::to_str(upper_boundaries[i], cvmodule->cv_width, cvmodule->cv_prec)+
402 ") is not commensurate to its bin width("+
403 cvm::to_str(widths[i], cvmodule->cv_width, cvmodule->cv_prec)+").\n");
404 upper_boundaries[i].real_value = lower_boundaries[i].real_value +
405 (nbins_round * widths[i]);
406 }
407
408 if (cvm::debug())
409 cvmodule->log("Number of points is "+cvm::to_str((int) nbins_round)+
410 " for the colvar no. "+cvm::to_str(i+1)+".\n");
411
412 nx.push_back(nbins_round);
413 }
414
415 return COLVARS_OK;
416 }
417
420 inline void wrap(std::vector<int> & ix) const
421 {
422 for (size_t i = 0; i < nd; i++) {
423 if (periodic[i]) {
424 ix[i] = (ix[i] + nx[i]) % nx[i]; // Avoid modulo with negative operands (implementation-defined)
425 } else {
426 if (ix[i] < 0 || ix[i] >= nx[i]) {
427 cvm::error_static("Trying to wrap illegal index vector (non-PBC) for a grid point: "
428 + cvm::to_str(ix), COLVARS_BUG_ERROR);
429 return;
430 }
431 }
432 }
433 }
434
437 inline bool wrap_detect_edge(std::vector<int> & ix) const
438 {
439 bool edge = false;
440 for (size_t i = 0; i < nd; i++) {
441 if (periodic[i]) {
442 ix[i] = (ix[i] + nx[i]) % nx[i]; // Avoid modulo with negative operands (implementation-defined)
443 } else if (ix[i] < 0 || ix[i] >= nx[i]) {
444 edge = true;
445 }
446 }
447 return edge;
448 }
449
452 inline bool wrap_to_edge(std::vector<int> & ix, std::vector<int> & edge_bin) const
453 {
454 bool edge = false;
455 edge_bin = ix;
456 for (size_t i = 0; i < nd; i++) {
457 if (periodic[i]) {
458 ix[i] = (ix[i] + nx[i]) % nx[i]; // Avoid modulo with negative operands (implementation-defined)
459 edge_bin[i] = ix[i];
460 } else if (ix[i] < 0) {
461 edge = true;
462 edge_bin[i] = 0;
463 } else if (ix[i] >= nx[i]) {
464 edge = true;
465 edge_bin[i] = nx[i] - 1;
466 }
467 }
468 return edge;
469 }
470
471
473 inline int current_bin_scalar(int const i) const
474 {
475 return value_to_bin_scalar(use_actual_value[i] ? cv[i]->actual_value() : cv[i]->value(), i);
476 }
477
480 inline int current_bin_flat_bound() const
481 {
482 std::vector<int> index = new_index();
483 for (size_t i = 0; i < nd; i++) {
484 index[i] = current_bin_scalar_bound(i);
485 }
486 return address(index);
487 }
488
491 inline int current_bin_scalar_bound(int const i) const
492 {
493 return value_to_bin_scalar_bound(use_actual_value[i] ? cv[i]->actual_value() : cv[i]->value(), i);
494 }
495
497 inline int current_bin_scalar(int const i, int const iv) const
498 {
500 cv[i]->actual_value().vector1d_value[iv] :
501 cv[i]->value().vector1d_value[iv], i);
502 }
503
506 inline int value_to_bin_scalar(colvarvalue const &value, const int i) const
507 {
508 return (int) cvm::floor( (value.real_value - lower_boundaries[i].real_value) / widths[i] );
509 }
510
512 inline cvm::real current_bin_scalar_fraction(int const i) const
513 {
514 return value_to_bin_scalar_fraction(use_actual_value[i] ? cv[i]->actual_value() : cv[i]->value(), i);
515 }
516
519 inline cvm::real value_to_bin_scalar_fraction(colvarvalue const &value, const int i) const
520 {
521 cvm::real x = (value.real_value - lower_boundaries[i].real_value) / widths[i];
522 return x - cvm::floor(x);
523 }
524
527 inline int value_to_bin_scalar_bound(colvarvalue const &value, const int i) const
528 {
529 int bin_index = cvm::floor( (value.real_value - lower_boundaries[i].real_value) / widths[i] );
530
531 // Wrap bins for periodic dimensions before truncating
532 if (periodic[i]) bin_index %= nx[i];
533 if (bin_index < 0) bin_index=0;
534 if (bin_index >=int(nx[i])) bin_index=int(nx[i])-1;
535 return (int) bin_index;
536 }
537
539 inline int value_to_bin_scalar(colvarvalue const &value,
540 colvarvalue const &new_offset,
541 cvm::real const &new_width) const
542 {
543 return (int) cvm::floor( (value.real_value - new_offset.real_value) / new_width );
544 }
545
548 inline colvarvalue bin_to_value_scalar(int const &i_bin, int const i) const
549 {
550 return lower_boundaries[i].real_value + widths[i] * (0.5 + i_bin);
551 }
552
554 inline colvarvalue bin_to_value_scalar(int const &i_bin,
555 colvarvalue const &new_offset,
556 cvm::real const &new_width) const
557 {
558 return new_offset.real_value + new_width * (0.5 + i_bin);
559 }
560
562 inline void set_value(std::vector<int> const &ix,
563 T const &t,
564 size_t const &imult = 0)
565 {
566 data[this->address(ix)+imult] = t;
567 has_data = true;
568 }
569
571 inline void set_value(size_t i, T const &t)
572 {
573 data[i] = t;
574 }
575
577 inline T get_value(size_t i) const
578 {
579 return data[i];
580 }
581
582
587 void delta_grid(colvar_grid<T> const &other_grid)
588 {
589
590 if (other_grid.multiplicity() != this->multiplicity()) {
591 cvm::error_static("Error: trying to subtract two grids with "
592 "different multiplicity.\n");
593 return;
594 }
595
596 if (other_grid.data.size() != this->data.size()) {
597 cvm::error_static("Error: trying to subtract two grids with "
598 "different size.\n");
599 return;
600 }
601
602 for (size_t i = 0; i < data.size(); i++) {
603 data[i] = other_grid.data[i] - data[i];
604 }
605 has_data = true;
606 }
607
608
612 void copy_grid(colvar_grid<T> const &other_grid)
613 {
614 if (other_grid.multiplicity() != this->multiplicity()) {
615 cvm::error_static("Error: trying to copy two grids with "
616 "different multiplicity.\n");
617 return;
618 }
619
620 if (other_grid.data.size() != this->data.size()) {
621 cvm::error_static("Error: trying to copy two grids with "
622 "different size.\n");
623 return;
624 }
625
626
627 for (size_t i = 0; i < data.size(); i++) {
628 data[i] = other_grid.data[i];
629 }
630 has_data = true;
631 }
632
635 void raw_data_out(T* out_data) const
636 {
637 for (size_t i = 0; i < data.size(); i++) out_data[i] = data[i];
638 }
639 void raw_data_out(std::vector<T>& out_data) const
640 {
641 out_data = data;
642 }
644 void raw_data_in(const T* in_data)
645 {
646 for (size_t i = 0; i < data.size(); i++) data[i] = in_data[i];
647 has_data = true;
648 }
649 void raw_data_in(const std::vector<T>& in_data)
650 {
651 data = in_data;
652 has_data = true;
653 }
655 size_t raw_data_num() const { return data.size(); }
656
657
660 inline T const & value(std::vector<int> const &ix,
661 size_t const &imult = 0) const
662 {
663 return data[this->address(ix) + imult];
664 }
665
667 inline T const & value(size_t i) const
668 {
669 return data[i];
670 }
671
673 inline void add_constant(T const &t)
674 {
675 for (size_t i = 0; i < nt; i++)
676 data[i] += t;
677 has_data = true;
678 }
679
681 inline void multiply_constant(cvm::real const &a)
682 {
683 for (size_t i = 0; i < nt; i++)
684 data[i] *= a;
685 }
686
688 inline void remove_small_values(cvm::real const &a)
689 {
690 for (size_t i = 0; i < nt; i++)
691 if(data[i]<a) data[i] = a;
692 }
693
694
697 inline std::vector<int> const get_colvars_index(std::vector<colvarvalue> const &values) const
698 {
699 std::vector<int> index = new_index();
700 for (size_t i = 0; i < nd; i++) {
701 index[i] = value_to_bin_scalar(values[i], i);
702 }
703 return index;
704 }
705
708 inline std::vector<int> const get_colvars_index() const
709 {
710 std::vector<int> index = new_index();
711 for (size_t i = 0; i < nd; i++) {
712 index[i] = current_bin_scalar(i);
713 }
714 return index;
715 }
716
719 inline std::vector<int> const get_colvars_index_bound() const
720 {
721 std::vector<int> index = new_index();
722 for (size_t i = 0; i < nd; i++) {
723 index[i] = current_bin_scalar_bound(i);
724 }
725 return index;
726 }
727
731 inline cvm::real bin_distance_from_boundaries(std::vector<colvarvalue> const &values,
732 bool skip_hard_boundaries = false)
733 {
734 cvm::real minimum = 1.0E+16;
735 for (size_t i = 0; i < nd; i++) {
736
737 if (periodic[i]) continue;
738
739 cvm::real dl = cvm::sqrt(cv[i]->dist2(values[i], lower_boundaries[i])) / widths[i];
740 cvm::real du = cvm::sqrt(cv[i]->dist2(values[i], upper_boundaries[i])) / widths[i];
741
742 if (values[i].real_value < lower_boundaries[i])
743 dl *= -1.0;
744 if (values[i].real_value > upper_boundaries[i])
745 du *= -1.0;
746
747 if ( ((!skip_hard_boundaries) || (!hard_lower_boundaries[i])) && (dl < minimum))
748 minimum = dl;
749 if ( ((!skip_hard_boundaries) || (!hard_upper_boundaries[i])) && (du < minimum))
750 minimum = du;
751 }
752
753 return minimum;
754 }
755
756
761 void map_grid(colvar_grid<T> const &other_grid)
762 {
763 if (other_grid.multiplicity() != this->multiplicity()) {
764 cvm::error_static("Error: trying to merge two grids with values of "
765 "different multiplicity.\n");
766 return;
767 }
768
769 std::vector<colvarvalue> const &gb = this->lower_boundaries;
770 std::vector<cvm::real> const &gw = this->widths;
771 std::vector<colvarvalue> const &ogb = other_grid.lower_boundaries;
772 std::vector<cvm::real> const &ogw = other_grid.widths;
773
774 std::vector<int> ix = this->new_index();
775 std::vector<int> oix = other_grid.new_index();
776
777 if (cvm::debug())
778 cvmodule->log("Remapping grid...\n");
779 for ( ; this->index_ok(ix); this->incr(ix)) {
780
781 for (size_t i = 0; i < nd; i++) {
782 oix[i] =
783 value_to_bin_scalar(bin_to_value_scalar(ix[i], gb[i], gw[i]),
784 ogb[i],
785 ogw[i]);
786 }
787
788 if (! other_grid.index_ok(oix)) {
789 continue;
790 }
791
792 for (size_t im = 0; im < mult; im++) {
793 this->set_value(ix, other_grid.value(oix, im), im);
794 }
795 }
796
797 has_data = true;
798 if (cvm::debug())
799 cvmodule->log("Remapping done.\n");
800 }
801
804 void add_grid(colvar_grid<T> const &other_grid,
805 cvm::real scale_factor = 1.0)
806 {
807 if (other_grid.multiplicity() != this->multiplicity()) {
808 cvm::error_static("Error: trying to sum togetehr two grids with values of "
809 "different multiplicity.\n");
810 return;
811 }
812 if (scale_factor != 1.0)
813 for (size_t i = 0; i < data.size(); i++) {
814 data[i] += static_cast<T>(scale_factor * other_grid.data[i]);
815 }
816 else
817 // skip multiplication if possible
818 for (size_t i = 0; i < data.size(); i++) {
819 data[i] += other_grid.data[i];
820 }
821 has_data = true;
822 }
823
826 virtual T value_output(std::vector<int> const &ix,
827 size_t const &imult = 0) const
828 {
829 return value(ix, imult);
830 }
831
835 virtual void value_input(std::vector<int> const &ix,
836 T const &t,
837 size_t const &imult = 0,
838 bool add = false)
839 {
840 if ( add )
841 data[address(ix) + imult] += t;
842 else
843 data[address(ix) + imult] = t;
844 has_data = true;
845 }
846
847
848 // /// Get the pointer to the binned value indexed by ix
849 // inline T const *value_p (std::vector<int> const &ix)
850 // {
851 // return &(data[address (ix)]);
852 // }
853
856 inline std::vector<int> const new_index() const
857 {
858 return std::vector<int> (nd, 0);
859 }
860
863 inline bool index_ok(std::vector<int> const &ix) const
864 {
865 for (size_t i = 0; i < nd; i++) {
866 if ( (ix[i] < 0) || (ix[i] >= int(nx[i])) )
867 return false;
868 }
869 return true;
870 }
871
874 inline void incr(std::vector<int> &ix) const
875 {
876 for (int i = ix.size()-1; i >= 0; i--) {
877
878 ix[i]++;
879
880 if (ix[i] >= nx[i]) {
881
882 if (i > 0) {
883 ix[i] = 0;
884 continue;
885 } else {
886 // this is the last iteration, a non-valid index is being
887 // set for the outer index, which will be caught by
888 // index_ok()
889 ix[0] = nx[0];
890 return;
891 }
892 } else {
893 return;
894 }
895 }
896 }
897
899 std::string get_state_params() const;
900
902 int parse_params(std::string const &conf,
904
909 {
910 for (size_t i = 0; i < nd; i++) {
911 if ( (cvm::sqrt(cv[i]->dist2(cv[i]->lower_boundary,
912 lower_boundaries[i])) > 1.0E-10) ||
913 (cvm::sqrt(cv[i]->dist2(cv[i]->upper_boundary,
914 upper_boundaries[i])) > 1.0E-10) ||
915 (cvm::sqrt(cv[i]->dist2(cv[i]->width,
916 widths[i])) > 1.0E-10) ) {
917 cvm::error_static("Error: restart information for a grid is "
918 "inconsistent with that of its colvars.\n");
919 return;
920 }
921 }
922 }
923
924
927 void check_consistency(colvar_grid<T> const &other_grid)
928 {
929 for (size_t i = 0; i < nd; i++) {
930 // we skip dist2(), because periodicities and the like should
931 // matter: boundaries should be EXACTLY the same (otherwise,
932 // map_grid() should be used)
933 if ( (cvm::fabs(other_grid.lower_boundaries[i] -
934 lower_boundaries[i]) > 1.0E-10) ||
935 (cvm::fabs(other_grid.upper_boundaries[i] -
936 upper_boundaries[i]) > 1.0E-10) ||
937 (cvm::fabs(other_grid.widths[i] -
938 widths[i]) > 1.0E-10) ||
939 (data.size() != other_grid.data.size()) ) {
940 cvm::error_static("Error: inconsistency between "
941 "two grids that are supposed to be equal, "
942 "aside from the data stored.\n");
943 return;
944 }
945 }
946 }
947
949 std::istream & read_restart(std::istream &is);
950
953
955 std::ostream & write_restart(std::ostream &os);
956
959
961 std::istream &read_raw(std::istream &is);
962
965
966private:
967
969 template <class IST> IST& read_restart_template_(IST& is);
970
972 template <class IST> IST& read_raw_template_(IST& is);
973
974public:
975
979 std::ostream &write_raw(std::ostream &os, size_t const buf_size = 3) const;
980
984 cvm::memory_stream &write_raw(cvm::memory_stream &os, size_t const buf_size = 3) const;
985
987 std::istream & read_multicol(std::istream &is, bool add = false);
988
990 int read_multicol(std::string const &filename,
991 std::string description = "grid file",
992 bool add = false);
993
995 std::ostream & write_multicol(std::ostream &os) const;
996
998 int write_multicol(std::string const &filename,
999 std::string description = "grid file") const;
1000
1002 std::ostream & write_opendx(std::ostream &os) const;
1003
1005 int write_opendx(std::string const &filename,
1006 std::string description = "grid file") const;
1007};
1008
1009
1010
1013class colvar_grid_count : public colvar_grid<size_t>
1014{
1015public:
1016
1019
1022 {}
1023
1025 colvar_grid_count(std::vector<colvar *> &colvars,
1026 std::shared_ptr<const colvar_grid_params> params = nullptr);
1027
1028 colvar_grid_count(std::vector<colvar *> &colvars,
1029 std::string config);
1030
1032 inline void incr_count(std::vector<int> const &ix)
1033 {
1034 ++(data[this->address(ix)]);
1035 }
1036
1038 inline size_t const & new_value(std::vector<int> const &ix)
1039 {
1040 return new_data[address(ix)];
1041 }
1042
1044 std::string get_state_params() const;
1045
1047 int parse_params(std::string const &conf,
1049
1051 std::istream & read_restart(std::istream &is);
1052
1055
1057 std::ostream & write_restart(std::ostream &os);
1058
1061
1063 std::istream &read_raw(std::istream &is);
1064
1067
1071 std::ostream &write_raw(std::ostream &os, size_t const buf_size = 3) const;
1072
1076 cvm::memory_stream &write_raw(cvm::memory_stream &os, size_t const buf_size = 3) const;
1077
1079 std::istream & read_multicol(std::istream &is, bool add = false);
1080
1082 int read_multicol(std::string const &filename,
1083 std::string description = "grid file",
1084 bool add = false);
1085
1087 std::ostream & write_multicol(std::ostream &os) const;
1088
1090 int write_multicol(std::string const &filename,
1091 std::string description = "grid file") const;
1092
1094 std::ostream & write_opendx(std::ostream &os) const;
1095
1097 int write_opendx(std::string const &filename,
1098 std::string description = "grid file") const;
1099
1101 virtual void value_input(std::vector<int> const &ix,
1102 size_t const &t,
1103 size_t const &imult = 0,
1104 bool add = false)
1105 {
1106 (void) imult;
1107 if (add) {
1108 data[address(ix)] += t;
1109 if (this->has_parent_data) {
1110 // save newly read data for inputting parent grid
1111 new_data[address(ix)] = t;
1112 }
1113 } else {
1114 data[address(ix)] = t;
1115 }
1116 has_data = true;
1117 }
1118
1121 inline int local_sample_count(int radius)
1122 {
1123 std::vector<int> ix0 = new_index();
1124 std::vector<int> ix = new_index();
1125
1126 for (size_t i = 0; i < nd; i++) {
1127 ix0[i] = current_bin_scalar_bound(i);
1128 }
1129 if (radius < 1) {
1130 // Simple case: no averaging
1131 if (index_ok(ix0))
1132 return value(ix0);
1133 else
1134 return 0;
1135 }
1136 size_t count = 0;
1137 size_t nbins = 0;
1138 int i, j, k;
1139 bool edge;
1140 ix = ix0;
1141 // Treat each dimension separately to simplify code
1142 switch (nd)
1143 {
1144 case 1:
1145 for (i = -radius; i <= radius; i++) {
1146 ix[0] = ix0[0] + i;
1147 edge = wrap_detect_edge(ix);
1148 if (!edge) {
1149 nbins++;
1150 count += value(ix);
1151 }
1152 }
1153 break;
1154 case 2:
1155 for (i = -radius; i <= radius; i++) {
1156 ix[0] = ix0[0] + i;
1157 for (j = -radius; j <= radius; j++) {
1158 ix[1] = ix0[1] + j;
1159 edge = wrap_detect_edge(ix);
1160 if (!edge) {
1161 nbins++;
1162 count += value(ix);
1163 }
1164 }
1165 }
1166 break;
1167 case 3:
1168 for (i = -radius; i <= radius; i++) {
1169 ix[0] = ix0[0] + i;
1170 for (j = -radius; j <= radius; j++) {
1171 ix[1] = ix0[1] + j;
1172 for (k = -radius; k <= radius; k++) {
1173 ix[2] = ix0[2] + k;
1174 edge = wrap_detect_edge(ix);
1175 if (!edge) {
1176 nbins++;
1177 count += value(ix);
1178 }
1179 }
1180 }
1181 }
1182 break;
1183 default:
1184 cvm::error_static("Error: local_sample_count is not implemented for grids of dimension > 3", COLVARS_NOT_IMPLEMENTED);
1185 break;
1186 }
1187
1188 if (nbins)
1189 // Integer division - an error on the order of 1 doesn't matter
1190 return count / nbins;
1191 else
1192 return 0.0;
1193 }
1194
1195
1199 inline cvm::real log_gradient_finite_diff(const std::vector<int> &ix0,
1200 int n = 0, int offset = 0)
1201 {
1202 cvm::real A0, A1, A2;
1203 std::vector<int> ix = ix0;
1204
1205 // TODO this can be rewritten more concisely with wrap_edge()
1206 if (periodic[n]) {
1207 ix[n]--; wrap(ix);
1208 A0 = value(ix) + offset;
1209 ix = ix0;
1210 ix[n]++; wrap(ix);
1211 A1 = value(ix) + offset;
1212 if (A0 * A1 == 0) {
1213 return 0.; // can't handle empty bins
1214 } else {
1215 return (cvmodule->logn(A1) - cvmodule->logn(A0))
1216 / (widths[n] * 2.);
1217 }
1218 } else if (ix[n] > 0 && ix[n] < nx[n]-1) { // not an edge
1219 ix[n]--;
1220 A0 = value(ix) + offset;
1221 ix = ix0;
1222 ix[n]++;
1223 A1 = value(ix) + offset;
1224 if (A0 * A1 == 0) {
1225 return 0.; // can't handle empty bins
1226 } else {
1227 return (cvmodule->logn(A1) - cvmodule->logn(A0))
1228 / (widths[n] * 2.);
1229 }
1230 } else {
1231 // edge: use 2nd order derivative
1232 int increment = (ix[n] == 0 ? 1 : -1);
1233 // move right from left edge, or the other way around
1234 A0 = value(ix) + offset;
1235 ix[n] += increment; A1 = value(ix) + offset;
1236 ix[n] += increment; A2 = value(ix) + offset;
1237 if (A0 * A1 * A2 == 0) {
1238 return 0.; // can't handle empty bins
1239 } else {
1240 return (-1.5 * cvmodule->logn(A0) + 2. * cvmodule->logn(A1)
1241 - 0.5 * cvmodule->logn(A2)) * increment / widths[n];
1242 }
1243 }
1244 }
1245
1246
1250 inline cvm::real gradient_finite_diff(const std::vector<int> &ix0,
1251 int n = 0)
1252 {
1253 cvm::real A0, A1, A2;
1254 std::vector<int> ix = ix0;
1255
1256 // FIXME this can be rewritten more concisely with wrap_edge()
1257 if (periodic[n]) {
1258 ix[n]--; wrap(ix);
1259 A0 = value(ix);
1260 ix = ix0;
1261 ix[n]++; wrap(ix);
1262 A1 = value(ix);
1263 if (A0 * A1 == 0) {
1264 return 0.; // can't handle empty bins
1265 } else {
1266 return (A1 - A0) / (widths[n] * 2.);
1267 }
1268 } else if (ix[n] > 0 && ix[n] < nx[n]-1) { // not an edge
1269 ix[n]--;
1270 A0 = value(ix);
1271 ix = ix0;
1272 ix[n]++;
1273 A1 = value(ix);
1274 if (A0 * A1 == 0) {
1275 return 0.; // can't handle empty bins
1276 } else {
1277 return (A1 - A0) / (widths[n] * 2.);
1278 }
1279 } else {
1280 // edge: use 2nd order derivative
1281 int increment = (ix[n] == 0 ? 1 : -1);
1282 // move right from left edge, or the other way around
1283 A0 = value(ix);
1284 ix[n] += increment; A1 = value(ix);
1285 ix[n] += increment; A2 = value(ix);
1286 return (-1.5 * A0 + 2. * A1
1287 - 0.5 * A2) * increment / widths[n];
1288 }
1289 }
1290};
1291
1292
1294class colvar_grid_scalar : public colvar_grid<cvm::real>
1295{
1296public:
1297
1301
1304
1307
1309 virtual ~colvar_grid_scalar();
1310
1312 colvar_grid_scalar(std::vector<colvar *> &colvars,
1313 std::shared_ptr<const colvar_grid_params> params = nullptr,
1314 bool add_extra_bin = false,
1315 std::string config = std::string());
1316
1318 colvar_grid_scalar(std::string const &filename);
1319
1321 inline void acc_value(std::vector<int> const &ix,
1322 cvm::real const &new_value,
1323 size_t const &imult = 0)
1324 {
1325 (void) imult;
1326 // only legal value of imult here is 0
1327 data[address(ix)] += new_value;
1328 if (samples)
1329 samples->incr_count(ix);
1330 has_data = true;
1331 }
1332
1334 std::string get_state_params() const;
1335
1337 int parse_params(std::string const &conf,
1339
1341 std::istream & read_restart(std::istream &is);
1342
1345
1347 std::ostream & write_restart(std::ostream &os);
1348
1351
1353 std::istream &read_raw(std::istream &is);
1354
1357
1361 std::ostream &write_raw(std::ostream &os, size_t const buf_size = 3) const;
1362
1366 cvm::memory_stream &write_raw(cvm::memory_stream &os, size_t const buf_size = 3) const;
1367
1369 std::istream & read_multicol(std::istream &is, bool add = false);
1370
1372 int read_multicol(std::string const &filename,
1373 std::string description = "grid file",
1374 bool add = false);
1375
1377 std::ostream & write_multicol(std::ostream &os) const;
1378
1380 int write_multicol(std::string const &filename,
1381 std::string description = "grid file") const;
1382
1384 std::ostream & write_opendx(std::ostream &os) const;
1385
1387 int write_opendx(std::string const &filename,
1388 std::string description = "grid file") const;
1389
1394 inline void vector_gradient_finite_diff( const std::vector<int> &ix0, std::vector<cvm::real> &grad)
1395 {
1396 cvm::real A0, A1;
1397 std::vector<int> ix;
1398 size_t i, j, k, n;
1399
1400 if (nd == 2) {
1401 for (n = 0; n < 2; n++) {
1402 ix = ix0;
1403 A0 = value(ix);
1404 ix[n]++; wrap(ix);
1405 A1 = value(ix);
1406 ix[1-n]++; wrap(ix);
1407 A1 += value(ix);
1408 ix[n]--; wrap(ix);
1409 A0 += value(ix);
1410 grad[n] = 0.5 * (A1 - A0) / widths[n];
1411 }
1412 } else if (nd == 3) {
1413
1414 cvm::real p[8]; // potential values within cube, indexed in binary (4 i + 2 j + k)
1415 ix = ix0;
1416 int index = 0;
1417 for (i = 0; i<2; i++) {
1418 ix[1] = ix0[1];
1419 for (j = 0; j<2; j++) {
1420 ix[2] = ix0[2];
1421 for (k = 0; k<2; k++) {
1422 wrap(ix);
1423 p[index++] = value(ix);
1424 ix[2]++;
1425 }
1426 ix[1]++;
1427 }
1428 ix[0]++;
1429 }
1430
1431 // The following would be easier to read using binary literals
1432 // 100 101 110 111 000 001 010 011
1433 grad[0] = 0.25 * ((p[4] + p[5] + p[6] + p[7]) - (p[0] + p[1] + p[2] + p[3])) / widths[0];
1434 // 010 011 110 111 000 001 100 101
1435 grad[1] = 0.25 * ((p[2] + p[3] + p[6] + p[7]) - (p[0] + p[1] + p[4] + p[5])) / widths[1];
1436 // 001 011 101 111 000 010 100 110
1437 grad[2] = 0.25 * ((p[1] + p[3] + p[5] + p[7]) - (p[0] + p[2] + p[4] + p[6])) / widths[2];
1438 } else {
1439 cvm::error_static("Finite differences available in dimension 2 and 3 only.");
1440 }
1441 }
1442
1443
1447 inline cvm::real log_gradient_finite_diff(const std::vector<int> &ix0,
1448 int n = 0, int offset = 0)
1449 {
1450 cvm::real A0, A1, A2;
1451 std::vector<int> ix = ix0;
1452
1453 // TODO this can be rewritten more concisely with wrap_edge()
1454 if (periodic[n]) {
1455 ix[n]--; wrap(ix);
1456 A0 = value(ix) + offset;
1457 ix = ix0;
1458 ix[n]++; wrap(ix);
1459 A1 = value(ix) + offset;
1460 if (A0 * A1 == 0) {
1461 return 0.; // can't handle empty bins
1462 } else {
1463 return (cvmodule->logn(A1) - cvmodule->logn(A0))
1464 / (widths[n] * 2.);
1465 }
1466 } else if (ix[n] > 0 && ix[n] < nx[n]-1) { // not an edge
1467 ix[n]--;
1468 A0 = value(ix) + offset;
1469 ix = ix0;
1470 ix[n]++;
1471 A1 = value(ix) + offset;
1472 if (A0 * A1 == 0) {
1473 return 0.; // can't handle empty bins
1474 } else {
1475 return (cvmodule->logn(A1) - cvmodule->logn(A0))
1476 / (widths[n] * 2.);
1477 }
1478 } else {
1479 // edge: use 2nd order derivative
1480 int increment = (ix[n] == 0 ? 1 : -1);
1481 // move right from left edge, or the other way around
1482 A0 = value(ix) + offset;
1483 ix[n] += increment; A1 = value(ix) + offset;
1484 ix[n] += increment; A2 = value(ix) + offset;
1485 if (A0 * A1 * A2 == 0) {
1486 return 0.; // can't handle empty bins
1487 } else {
1488 return (-1.5 * cvmodule->logn(A0) + 2. * cvmodule->logn(A1)
1489 - 0.5 * cvmodule->logn(A2)) * increment / widths[n];
1490 }
1491 }
1492 }
1493
1494
1498 inline cvm::real gradient_finite_diff(const std::vector<int> &ix0,
1499 int n = 0)
1500 {
1501 cvm::real A0, A1, A2;
1502 std::vector<int> ix = ix0;
1503
1504 // FIXME this can be rewritten more concisely with wrap_edge()
1505 if (periodic[n]) {
1506 ix[n]--; wrap(ix);
1507 A0 = value(ix);
1508 ix = ix0;
1509 ix[n]++; wrap(ix);
1510 A1 = value(ix);
1511 if (A0 * A1 == 0) {
1512 return 0.; // can't handle empty bins
1513 } else {
1514 return (A1 - A0) / (widths[n] * 2.);
1515 }
1516 } else if (ix[n] > 0 && ix[n] < nx[n]-1) { // not an edge
1517 ix[n]--;
1518 A0 = value(ix);
1519 ix = ix0;
1520 ix[n]++;
1521 A1 = value(ix);
1522 if (A0 * A1 == 0) {
1523 return 0.; // can't handle empty bins
1524 } else {
1525 return cvm::real(A1 - A0) / (widths[n] * 2.);
1526 }
1527 } else {
1528 // edge: use 2nd order derivative
1529 int increment = (ix[n] == 0 ? 1 : -1);
1530 // move right from left edge, or the other way around
1531 A0 = value(ix);
1532 ix[n] += increment; A1 = value(ix);
1533 ix[n] += increment; A2 = value(ix);
1534 return (-1.5 * cvm::real(A0) + 2. * cvm::real(A1)
1535 - 0.5 * cvm::real(A2)) * increment / widths[n];
1536 }
1537 }
1538
1539
1542 virtual inline cvm::real value_output(std::vector<int> const &ix,
1543 size_t const &imult = 0) const override
1544 {
1545 int s;
1546 if (imult > 0) {
1547 cvm::error_static("Error: trying to access a component "
1548 "larger than 1 in a scalar data grid.\n");
1549 return 0.;
1550 }
1551 if (samples) {
1552 return ( (s = samples->value(ix)) > 0) ?
1553 (data[address(ix) + imult] / cvm::real(s)) :
1554 0.0;
1555 } else {
1556 return data[address(ix) + imult];
1557 }
1558 }
1559
1561 virtual void value_input(std::vector<int> const &ix,
1562 cvm::real const &new_value,
1563 size_t const &imult = 0,
1564 bool add = false) override
1565 {
1566 if (imult > 0) {
1567 cvm::error_static("Error: trying to access a component "
1568 "larger than 1 in a scalar data grid.\n");
1569 return;
1570 }
1571 if (add) {
1572 if (samples)
1573 data[address(ix)] += new_value * samples->new_value(ix);
1574 else
1575 data[address(ix)] += new_value;
1576 } else {
1577 if (samples)
1578 data[address(ix)] = new_value * samples->value(ix);
1579 else
1580 data[address(ix)] = new_value;
1581 }
1582 has_data = true;
1583 }
1584
1586 cvm::real maximum_value() const;
1587
1589 cvm::real minimum_value() const;
1590
1593
1595 cvm::real integral() const;
1596
1599 cvm::real entropy() const;
1600
1603 cvm::real grid_rmsd(colvar_grid_scalar const &other_grid) const;
1604};
1605
1606
1607
1609class colvar_grid_gradient : public colvar_grid<cvm::real>
1610{
1611public:
1612
1615 std::shared_ptr<colvar_grid_count> samples;
1616
1619
1622 {}
1623
1624 // /// Constructor from specific sizes arrays
1625 // colvar_grid_gradient(std::vector<int> const &nx_i);
1626
1627 // /// Constructor from a vector of colvars
1628 // colvar_grid_gradient(std::vector<colvar *> &colvars,
1629 // std::string config = std::string());
1630
1632 colvar_grid_gradient(std::string const &filename);
1633
1635 colvar_grid_gradient(std::vector<colvar *> &colvars,
1636 std::shared_ptr<colvar_grid_count> samples_in = nullptr,
1637 std::shared_ptr<const colvar_grid_params> params = nullptr,
1638 std::string config = std::string());
1639
1642 int min_samples;
1643
1645 std::string get_state_params() const;
1646
1648 int parse_params(std::string const &conf,
1650
1652 std::istream & read_restart(std::istream &is);
1653
1656
1658 std::ostream & write_restart(std::ostream &os);
1659
1662
1664 std::istream &read_raw(std::istream &is);
1665
1668
1672 std::ostream &write_raw(std::ostream &os, size_t const buf_size = 3) const;
1673
1677 cvm::memory_stream &write_raw(cvm::memory_stream &os, size_t const buf_size = 3) const;
1678
1680 std::istream & read_multicol(std::istream &is, bool add = false);
1681
1683 int read_multicol(std::string const &filename,
1684 std::string description = "grid file",
1685 bool add = false);
1686
1688 std::ostream & write_multicol(std::ostream &os) const;
1689
1691 int write_multicol(std::string const &filename,
1692 std::string description = "grid file") const;
1693
1695 std::ostream & write_opendx(std::ostream &os) const;
1696
1698 int write_opendx(std::string const &filename,
1699 std::string description = "grid file") const;
1700
1702 inline void vector_value(std::vector<int> const &ix, std::vector<cvm::real> &v) const
1703 {
1704 cvm::real const * p = &value(ix);
1705 if (samples) {
1706 int count = samples->value(ix);
1707 if (count) {
1708 cvm::real invcount = 1.0 / count;
1709 for (size_t i = 0; i < mult; i++) {
1710 v[i] = invcount * p[i];
1711 }
1712 } else {
1713 for (size_t i = 0; i < mult; i++) {
1714 v[i] = 0.0;
1715 }
1716 }
1717 } else {
1718 for (size_t i = 0; i < mult; i++) {
1719 v[i] = p[i];
1720 }
1721 }
1722 }
1723
1724
1726 inline void acc_value(std::vector<int> const &ix, std::vector<colvarvalue> const &values) {
1727 for (size_t imult = 0; imult < mult; imult++) {
1728 data[address(ix) + imult] += values[imult].real_value;
1729 }
1730 if (samples)
1731 samples->incr_count(ix);
1732 }
1733
1736 inline void acc_force(std::vector<int> const &ix, cvm::real const *forces) {
1737 for (size_t imult = 0; imult < mult; imult++) {
1738 data[address(ix) + imult] -= forces[imult];
1739 }
1740 if (samples)
1741 samples->incr_count(ix);
1742 }
1743
1746 virtual cvm::real value_output(std::vector<int> const &ix,
1747 size_t const &imult = 0) const override
1748 {
1749 int s;
1750 if (samples) {
1751 return ( (s = samples->value(ix)) > 0) ?
1752 (data[address(ix) + imult] / cvm::real(s)) :
1753 0.0;
1754 } else {
1755 return data[address(ix) + imult];
1756 }
1757 }
1758
1762 {
1763 cvm::real fact;
1764 if ( weight <= min_samples ) {
1765 fact = 0.0;
1766 } else if ( weight < full_samples ) {
1767 fact = (weight - min_samples) / (weight * cvm::real(full_samples - min_samples));
1768 } else {
1769 fact = 1.0 / weight;
1770 }
1771 return fact;
1772 }
1773
1774
1779 virtual inline cvm::real value_output_smoothed(std::vector<int> const &ix, bool smoothed = true)
1780 {
1781 cvm::real weight, fact;
1782
1783 if (samples) {
1784 weight = cvm::real(samples->value(ix));
1785 } else {
1786 weight = 1.;
1787 }
1788
1789 if (smoothed) {
1790 fact = smooth_inverse_weight(weight);
1791 } else {
1792 fact = weight > 0. ? 1. / weight : 0.;
1793 }
1794
1795 return fact * data[address(ix)];
1796 }
1797
1801 inline void vector_value_smoothed(std::vector<int> const &ix, cvm::real *grad, bool smoothed = true)
1802 {
1803 cvm::real weight, fact;
1804
1805 if (samples) {
1806 weight = cvm::real(samples->value(ix));
1807 } else {
1808 weight = 1.;
1809 }
1810
1811 if (smoothed) {
1812 fact = smooth_inverse_weight(weight);
1813 } else {
1814 fact = weight > 0. ? 1. / weight : 0.;
1815 }
1816
1817 cvm::real *p = &(data[address(ix)]);
1818
1819 for (size_t imult = 0; imult < mult; imult++) {
1820 grad[imult] = fact * p[imult];
1821 }
1822 }
1823
1827 virtual void value_input(std::vector<int> const &ix,
1828 cvm::real const &new_value,
1829 size_t const &imult = 0,
1830 bool add = false) override
1831 {
1832 if (add) {
1833 if (samples)
1834 data[address(ix) + imult] += new_value * samples->new_value(ix);
1835 else
1836 data[address(ix) + imult] += new_value;
1837 } else {
1838 if (samples)
1839 data[address(ix) + imult] = new_value * samples->value(ix);
1840 else
1841 data[address(ix) + imult] = new_value;
1842 }
1843 has_data = true;
1844 }
1845
1846
1848 inline cvm::real average(bool smoothed = false)
1849 {
1850 if (nd != 1 || nx[0] == 0) {
1851 return 0.0;
1852 }
1853
1854 cvm::real sum = 0.0;
1855 for (std::vector<int> ix = new_index(); index_ok(ix); incr(ix)) {
1856 sum += value_output_smoothed(ix, smoothed);
1857 }
1858
1859 return (sum / cvm::real(nx[0]));
1860 }
1861
1864 cvm::real grid_rmsd(colvar_grid_gradient const &other_grid) const;
1865
1868 void write_1D_integral(std::ostream &os);
1869
1870};
1871
1872#endif
1873
Colvar_grid derived class to hold counters in discrete n-dim colvar space.
Definition: colvargrid.h:1014
colvar_grid_count()
Default constructor.
Definition: colvargrid.cpp:21
int local_sample_count(int radius)
Return the average number of samples in a given "radius" around current bin Really a hypercube of len...
Definition: colvargrid.h:1121
std::istream & read_raw(std::istream &is)
Read all grid parameters and data from a formatted stream.
Definition: colvargrid.cpp:68
std::string get_state_params() const
Write the current grid parameters to a string.
Definition: colvargrid.cpp:37
cvm::real gradient_finite_diff(const std::vector< int > &ix0, int n=0)
Return the gradient of discrete count from finite differences on the same grid for dimension n (colva...
Definition: colvargrid.h:1250
std::ostream & write_restart(std::ostream &os)
Write all grid parameters and data to a formatted stream.
Definition: colvargrid.cpp:58
std::ostream & write_opendx(std::ostream &os) const
Write the grid data without labels, as they are represented in memory.
Definition: colvargrid.cpp:112
std::ostream & write_raw(std::ostream &os, size_t const buf_size=3) const
Definition: colvargrid.cpp:78
int parse_params(std::string const &conf, colvarparse::Parse_Mode const parse_mode=colvarparse::parse_normal)
Read new grid parameters from a string.
Definition: colvargrid.cpp:42
virtual void value_input(std::vector< int > const &ix, size_t const &t, size_t const &imult=0, bool add=false)
Enter or add a value, but also handle parent grid.
Definition: colvargrid.h:1101
std::istream & read_restart(std::istream &is)
Read all grid parameters and data from a formatted stream.
Definition: colvargrid.cpp:48
std::istream & read_multicol(std::istream &is, bool add=false)
Read a grid written by write_multicol(), incrementin if data is true.
Definition: colvargrid.cpp:89
std::ostream & write_multicol(std::ostream &os) const
Write grid in a format which is both human-readable and gnuplot-friendly.
Definition: colvargrid.cpp:101
cvm::real log_gradient_finite_diff(const std::vector< int > &ix0, int n=0, int offset=0)
Return the log-gradient from finite differences on the same grid for dimension n (colvar_grid_count)
Definition: colvargrid.h:1199
size_t const & new_value(std::vector< int > const &ix)
Get the binned count indexed by ix from the newly read data.
Definition: colvargrid.h:1038
virtual ~colvar_grid_count()
Destructor.
Definition: colvargrid.h:1021
void incr_count(std::vector< int > const &ix)
Increment the counter at given position.
Definition: colvargrid.h:1032
Class for accumulating the gradient of a scalar function on a grid.
Definition: colvargrid.h:1610
std::string get_state_params() const
Write the current grid parameters to a string.
Definition: colvargrid.cpp:369
cvm::real smooth_inverse_weight(cvm::real weight)
Definition: colvargrid.h:1761
void vector_value_smoothed(std::vector< int > const &ix, cvm::real *grad, bool smoothed=true)
Obtain the vector value of the function at ix divided by its number of samples (if the count grid is ...
Definition: colvargrid.h:1801
virtual cvm::real value_output_smoothed(std::vector< int > const &ix, bool smoothed=true)
Return the scalar value of the function at ix divided by its number of samples (if the count grid is ...
Definition: colvargrid.h:1779
void acc_value(std::vector< int > const &ix, std::vector< colvarvalue > const &values)
Accumulate the value.
Definition: colvargrid.h:1726
cvm::real average(bool smoothed=false)
Compute and return average value for a 1D gradient grid.
Definition: colvargrid.h:1848
std::shared_ptr< colvar_grid_count > samples
Provide the sample count by which each binned value should be divided.
Definition: colvargrid.h:1615
std::istream & read_multicol(std::istream &is, bool add=false)
Read a grid written by write_multicol(), incrementin if data is true.
Definition: colvargrid.cpp:421
colvar_grid_gradient()
Default constructor.
Definition: colvargrid.cpp:335
std::ostream & write_restart(std::ostream &os)
Write all grid parameters and data to a formatted stream.
Definition: colvargrid.cpp:390
void write_1D_integral(std::ostream &os)
If the grid is 1-dimensional, integrate it and write the integral to a file (DEPRECATED by the colvar...
Definition: colvargrid.cpp:456
int parse_params(std::string const &conf, colvarparse::Parse_Mode const parse_mode=colvarparse::parse_normal)
Read new grid parameters from a string.
Definition: colvargrid.cpp:374
virtual cvm::real value_output(std::vector< int > const &ix, size_t const &imult=0) const override
Return the value of the function at ix divided by its number of samples (if the count grid is defined...
Definition: colvargrid.h:1746
std::istream & read_raw(std::istream &is)
Read all grid parameters and data from a formatted stream.
Definition: colvargrid.cpp:400
virtual ~colvar_grid_gradient()
Destructor.
Definition: colvargrid.h:1621
virtual void value_input(std::vector< int > const &ix, cvm::real const &new_value, size_t const &imult=0, bool add=false) override
Get the value from a formatted output and transform it into the internal representation (it may have ...
Definition: colvargrid.h:1827
void acc_force(std::vector< int > const &ix, cvm::real const *forces)
Accumulate the gradient based on the force (i.e. sums the opposite of the force)
Definition: colvargrid.h:1736
int full_samples
Parameters for smoothing data with low sampling.
Definition: colvargrid.h:1641
std::ostream & write_raw(std::ostream &os, size_t const buf_size=3) const
Definition: colvargrid.cpp:410
std::ostream & write_multicol(std::ostream &os) const
Write grid in a format which is both human-readable and gnuplot-friendly.
Definition: colvargrid.cpp:433
std::istream & read_restart(std::istream &is)
Read all grid parameters and data from a formatted stream.
Definition: colvargrid.cpp:380
cvm::real grid_rmsd(colvar_grid_gradient const &other_grid) const
Return the RMSD between this grid's data and another one Grids must have the same dimensions.
Definition: colvargrid.cpp:513
void vector_value(std::vector< int > const &ix, std::vector< cvm::real > &v) const
Get a vector with the binned value(s) indexed by ix, normalized if applicable.
Definition: colvargrid.h:1702
std::ostream & write_opendx(std::ostream &os) const
Write the grid data without labels, as they are represented in memory.
Definition: colvargrid.cpp:444
Unified base class for grid of values of a function of several collective variables.
Definition: colvargrid.h:24
size_t nd
Number of dimensions.
Definition: colvargrid.h:28
std::vector< int > nxc
Cumulative number of points along each dimension.
Definition: colvargrid.h:34
std::vector< colvarvalue > lower_boundaries
Lower boundaries of the colvars in this grid.
Definition: colvargrid.h:37
std::vector< cvm::real > widths
Widths of the colvars in this grid.
Definition: colvargrid.h:43
std::vector< colvarvalue > upper_boundaries
Upper boundaries of the colvars in this grid.
Definition: colvargrid.h:40
std::vector< int > nx
Number of points along each dimension.
Definition: colvargrid.h:31
Class for accumulating a scalar function on a grid.
Definition: colvargrid.h:1295
void acc_value(std::vector< int > const &ix, cvm::real const &new_value, size_t const &imult=0)
Accumulate the value.
Definition: colvargrid.h:1321
cvm::real minimum_value() const
Return the lowest value.
Definition: colvargrid.cpp:249
cvm::real gradient_finite_diff(const std::vector< int > &ix0, int n=0)
Return the gradient of discrete count from finite differences on the same grid for dimension n (colva...
Definition: colvargrid.h:1498
std::istream & read_restart(std::istream &is)
Read all grid parameters and data from a formatted stream.
Definition: colvargrid.cpp:163
std::string get_state_params() const
Write the current grid parameters to a string.
Definition: colvargrid.cpp:152
virtual ~colvar_grid_scalar()
Destructor.
Definition: colvargrid.cpp:148
cvm::real minimum_pos_value() const
Return the lowest positive value.
Definition: colvargrid.cpp:258
virtual cvm::real value_output(std::vector< int > const &ix, size_t const &imult=0) const override
Return the value of the function at ix divided by its number of samples (if the count grid is defined...
Definition: colvargrid.h:1542
cvm::real grid_rmsd(colvar_grid_scalar const &other_grid) const
Return the RMSD between this grid's data and another one Grids must have the same dimensions.
Definition: colvargrid.cpp:305
std::ostream & write_opendx(std::ostream &os) const
Write the grid data without labels, as they are represented in memory.
Definition: colvargrid.cpp:227
std::ostream & write_restart(std::ostream &os)
Write all grid parameters and data to a formatted stream.
Definition: colvargrid.cpp:173
std::ostream & write_multicol(std::ostream &os) const
Write grid in a format which is both human-readable and gnuplot-friendly.
Definition: colvargrid.cpp:216
std::istream & read_multicol(std::istream &is, bool add=false)
Read a grid written by write_multicol(), incrementin if data is true.
Definition: colvargrid.cpp:204
cvm::real integral() const
Calculates the integral of the map (uses widths if they are defined)
Definition: colvargrid.cpp:274
int parse_params(std::string const &conf, colvarparse::Parse_Mode const parse_mode=colvarparse::parse_normal)
Read new grid parameters from a string.
Definition: colvargrid.cpp:157
std::ostream & write_raw(std::ostream &os, size_t const buf_size=3) const
Definition: colvargrid.cpp:193
virtual void value_input(std::vector< int > const &ix, cvm::real const &new_value, size_t const &imult=0, bool add=false) override
Enter or add value but also deal with count grid.
Definition: colvargrid.h:1561
cvm::real entropy() const
Assuming that the map is a normalized probability density, calculates the entropy (uses widths if the...
Definition: colvargrid.cpp:288
colvar_grid_scalar()
Default constructor.
Definition: colvargrid.cpp:125
cvm::real maximum_value() const
Return the highest value.
Definition: colvargrid.cpp:239
cvm::real log_gradient_finite_diff(const std::vector< int > &ix0, int n=0, int offset=0)
Return the log-gradient from finite differences on the same grid for dimension n (colvar_grid_scalar)
Definition: colvargrid.h:1447
colvar_grid_count * samples
Provide the associated sample count by which each binned value should be divided.
Definition: colvargrid.h:1300
void vector_gradient_finite_diff(const std::vector< int > &ix0, std::vector< cvm::real > &grad)
Return the gradient of the scalar field from finite differences Input coordinates are those of gradie...
Definition: colvargrid.h:1394
std::istream & read_raw(std::istream &is)
Read all grid parameters and data from a formatted stream.
Definition: colvargrid.cpp:183
Grid of values of a function of several collective variables.
Definition: colvargrid.h:52
int setup()
Allocate data (allow initialization also after construction)
Definition: colvargrid.h:200
std::istream & read_restart(std::istream &is)
Read all grid parameters and data from a formatted stream.
Definition: colvargrid_def.h:102
void set_value(size_t i, T const &t)
Set the value at the point with linear address i (for speed)
Definition: colvargrid.h:571
void check_consistency()
Check that the grid information inside (boundaries, widths, ...) is consistent with the current setti...
Definition: colvargrid.h:908
std::string get_state_params() const
Write the current grid parameters to a string.
Definition: colvargrid_def.h:203
std::ostream & write_raw(std::ostream &os, size_t const buf_size=3) const
Definition: colvargrid_def.h:170
int current_bin_scalar_bound(int const i) const
Report the bin corresponding to the current value of variable i and assign first or last bin if out o...
Definition: colvargrid.h:491
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 loopi...
Definition: colvargrid.h:856
void add_constant(T const &t)
Add a constant to all elements (fast loop)
Definition: colvargrid.h:673
std::vector< size_t > new_data
Newly read data (used for count grids, when adding several grids read from disk)
Definition: colvargrid.h:68
cvm::real value_to_bin_scalar_fraction(colvarvalue const &value, const int i) const
Use the lower boundary and the width to report the fraction of bin beyond value_to_bin_scalar() that ...
Definition: colvargrid.h:519
void incr(std::vector< int > &ix) const
Increment the index, in a way that will make it loop over the whole nd-dimensional array.
Definition: colvargrid.h:874
size_t number_of_points(int const icv=-1) const
Definition: colvargrid.h:122
std::vector< bool > hard_upper_boundaries
Whether some colvars have hard upper boundaries.
Definition: colvargrid.h:100
virtual void value_input(std::vector< int > const &ix, T const &t, size_t const &imult=0, bool add=false)
Get the value from a formatted output and transform it into the internal representation (the two may ...
Definition: colvargrid.h:835
T const & value(size_t i) const
Get the binned value indexed by linear address i.
Definition: colvargrid.h:667
std::vector< int > const get_colvars_index_bound() const
Get the bin indices corresponding to the provided values of the colvars and assign first or last bin ...
Definition: colvargrid.h:719
int current_bin_flat_bound() const
Report the flattened bin address corresponding to the current value of all variables and assign first...
Definition: colvargrid.h:480
size_t mult
Multiplicity of each datum (allow the binning of non-scalar types such as atomic gradients)
Definition: colvargrid.h:59
size_t multiplicity() const
Return the multiplicity of the type used.
Definition: colvargrid.h:144
std::vector< int > const & sizes() const
Get the sizes in each direction.
Definition: colvargrid.h:132
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.
Definition: colvargrid.h:562
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.
Definition: colvargrid.h:548
std::vector< T > data
Low-level array of values.
Definition: colvargrid.h:65
void multiply_constant(cvm::real const &a)
Multiply all elements by a scalar constant (fast loop)
Definition: colvargrid.h:681
std::vector< int > const & number_of_points_vec() const
Return the numbers of points in all dimensions.
Definition: colvargrid.h:115
bool index_ok(std::vector< int > const &ix) const
Check that the index is within range in each of the dimensions.
Definition: colvargrid.h:863
cvm::real current_bin_scalar_fraction(int const i) const
Report the fraction of bin beyond current_bin_scalar()
Definition: colvargrid.h:512
bool has_data
Whether this grid has been filled with data or is still empty.
Definition: colvargrid.h:106
void add_grid(colvar_grid< T > const &other_grid, cvm::real scale_factor=1.0)
Add data from another grid of the same type, AND identical definition (boundaries,...
Definition: colvargrid.h:804
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.
Definition: colvargrid.h:697
void check_consistency(colvar_grid< T > const &other_grid)
Check that the grid information inside (boundaries, widths, ...) is consistent with that of another g...
Definition: colvargrid.h:927
colvar_grid(std::vector< int > const &nx_i, T const &t=T(), size_t mult_i=1)
Constructor from explicit grid sizes.
Definition: colvargrid.h:247
int parse_params(std::string const &conf, colvarparse::Parse_Mode const parse_mode=colvarparse::parse_normal)
Read new grid parameters from a string.
Definition: colvargrid_def.h:233
size_t nt
Total number of grid points.
Definition: colvargrid.h:62
T get_value(size_t i) const
Get the value at the point with linear address i (for speed)
Definition: colvargrid.h:577
void reset(T const &t=T())
Reset data (in case the grid is being reused)
Definition: colvargrid.h:206
size_t num_variables() const
Return the number of colvar objects.
Definition: colvargrid.h:109
int setup(std::vector< int > const &nx_i, T const &t=T(), size_t const &mult_i=1)
Allocate data.
Definition: colvargrid.h:159
IST & read_restart_template_(IST &is)
Helper method for read_restart.
Definition: colvargrid_def.h:84
colvar_grid(std::vector< colvar * > const &colvars, T const &t=T(), size_t mult_i=1, bool add_extra_bin=false, std::shared_ptr< const colvar_grid_params > params=nullptr, std::string config=std::string())
Constructor from a vector of colvars or an optional grid config string.
Definition: colvargrid.h:258
colvar_grid(colvar_grid< T > const &g)
"Almost copy-constructor": only copies configuration parameters from another grid,...
Definition: colvargrid.h:230
bool has_parent_data
True if this is a count grid related to another grid of data.
Definition: colvargrid.h:103
std::ostream & write_opendx(std::ostream &os) const
Write the grid data without labels, as they are represented in memory.
Definition: colvargrid_def.h:494
bool wrap_to_edge(std::vector< int > &ix, std::vector< int > &edge_bin) const
Definition: colvargrid.h:452
std::vector< int > const get_colvars_index() const
Get the bin indices corresponding to the current values of the colvars.
Definition: colvargrid.h:708
int current_bin_scalar(int const i, int const iv) const
Report the bin corresponding to the current value of item iv in variable i.
Definition: colvargrid.h:497
void delta_grid(colvar_grid< T > const &other_grid)
Get the change from this to other_grid and store the result in this. this_grid := other_grid - this_g...
Definition: colvargrid.h:587
void request_actual_value(bool b=true)
Request grid to use actual values of extended coords.
Definition: colvargrid.h:150
cvm::real bin_distance_from_boundaries(std::vector< colvarvalue > const &values, bool skip_hard_boundaries=false)
Get the minimal distance (in number of bins) from the boundaries; a negative number is returned if th...
Definition: colvargrid.h:731
int value_to_bin_scalar_bound(colvarvalue const &value, const int i) const
Use the lower boundary and the width to report which bin the provided value is in and assign first or...
Definition: colvargrid.h:527
int current_bin_scalar(int const i) const
Report the bin corresponding to the current value of variable i.
Definition: colvargrid.h:473
std::istream & read_raw(std::istream &is)
Read all grid parameters and data from a formatted stream.
Definition: colvargrid_def.h:157
std::vector< bool > hard_lower_boundaries
Whether some colvars have hard lower boundaries.
Definition: colvargrid.h:97
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.
Definition: colvargrid.h:660
std::vector< bool > use_actual_value
Do we request actual value (for extended-system colvars)?
Definition: colvargrid.h:74
size_t raw_data_num() const
Size of the data as they are represented in memory.
Definition: colvargrid.h:655
size_t address(std::vector< int > const &ix) const
Get the low-level index corresponding to an index.
Definition: colvargrid.h:77
void remove_small_values(cvm::real const &a)
Assign values that are smaller than scalar constant the latter value (fast loop)
Definition: colvargrid.h:688
std::ostream & write_multicol(std::ostream &os) const
Write grid in a format which is both human-readable and gnuplot-friendly.
Definition: colvargrid_def.h:430
virtual ~colvar_grid()
Destructor.
Definition: colvargrid.h:223
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.
Definition: colvargrid.h:506
IST & read_raw_template_(IST &is)
Helper method for read_raw.
Definition: colvargrid_def.h:129
virtual T value_output(std::vector< int > const &ix, size_t const &imult=0) const
Return the value suitable for output purposes (so that it may be rescaled or manipulated without chan...
Definition: colvargrid.h:826
void set_sizes(std::vector< int > const &new_sizes)
Set the sizes in each direction.
Definition: colvargrid.h:138
void copy_grid(colvar_grid< T > const &other_grid)
Copy data from another grid of the same type, AND identical definition (boundaries,...
Definition: colvargrid.h:612
bool wrap_detect_edge(std::vector< int > &ix) const
Definition: colvargrid.h:437
void raw_data_out(T *out_data) const
Extract the grid data as they are represented in memory. Put the results in "out_data".
Definition: colvargrid.h:635
void map_grid(colvar_grid< T > const &other_grid)
Add data from another grid of the same type.
Definition: colvargrid.h:761
std::istream & read_multicol(std::istream &is, bool add=false)
Read a grid written by write_multicol(), incrementing if add is true.
Definition: colvargrid_def.h:311
std::vector< bool > periodic
Whether some colvars are periodic.
Definition: colvargrid.h:94
int value_to_bin_scalar(colvarvalue const &value, colvarvalue const &new_offset, cvm::real const &new_width) const
Same as the standard version, but uses another grid definition.
Definition: colvargrid.h:539
void wrap(std::vector< int > &ix) const
Definition: colvargrid.h:420
colvarvalue bin_to_value_scalar(int const &i_bin, colvarvalue const &new_offset, cvm::real const &new_width) const
Same as the standard version, but uses different parameters.
Definition: colvargrid.h:554
colvar_grid()
Default constructor.
Definition: colvargrid.h:214
std::vector< colvar * > cv
Colvars collected in this grid.
Definition: colvargrid.h:71
std::ostream & write_restart(std::ostream &os)
Write all grid parameters and data to a formatted stream.
Definition: colvargrid_def.h:112
void raw_data_in(const T *in_data)
Input the data as they are represented in memory.
Definition: colvargrid.h:644
@ f_cv_hard_lower_boundary
The lower boundary is not defined from user's choice.
Definition: colvardeps.h:327
@ f_cv_hard_upper_boundary
The upper boundary is not defined from user's choice.
Definition: colvardeps.h:329
Collective variables module (main class)
Definition: colvarmodule.h:72
static real logn(real const &x)
Definition: colvarmodule.h:220
double real
Defining an abstract real number allows to switch precision.
Definition: colvarmodule.h:99
static real floor(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:125
static int error_static(std::string const &message, int code=-1)
Definition: colvarmodule.h:770
void log(std::string const &message, int min_log_level=10)
Definition: colvarmodule.cpp:2102
int error(std::string const &message, int code=-1)
Print a message to the main log and set global error code.
Definition: colvarmodule.cpp:2186
static constexpr size_t cv_prec
Number of digits to represent a collective variables value(s)
Definition: colvarmodule.h:724
static real sqrt(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:143
static bool debug()
Whether debug output should be enabled (compile-time option)
Definition: colvarmodule.h:376
static real fabs(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:137
static constexpr size_t cv_width
Number of characters to represent a collective variables value(s)
Definition: colvarmodule.h:726
static std::string to_str(char const *s)
Convert to string for output purposes.
Definition: colvarmodule.cpp:2541
Base class containing parsing functions; all objects which need to parse input inherit from this.
Definition: colvarparse.h:27
int check_keywords(std::string &conf, char const *key)
Check that all the keywords within "conf" are in the list of allowed keywords; this will invoke strip...
Definition: colvarparse.cpp:598
Parse_Mode
How a keyword is parsed in a string.
Definition: colvarparse.h:63
@ parse_normal
Alias for old default behavior (should be phased out)
Definition: colvarparse.h:82
Value of a collective variable: this is a metatype which can be set at runtime. By default it is set ...
Definition: colvarvalue.h:43
@ type_scalar
Scalar number, implemented as colvarmodule::real (default)
Definition: colvarvalue.h:56
cvm::real real_value
Real data member.
Definition: colvarvalue.h:77
Definition: colvars_memstream.h:30
Collective variables main module.
Parsing functions for collective variables.