Collective Variables Module - Developer Documentation
Loading...
Searching...
No Matches
colvar.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 COLVAR_H
11#define COLVAR_H
12
13#include <iostream>
14#include <map>
15#include <functional>
16
17#include "colvarmodule.h"
18#include "colvarvalue.h"
19#include "colvarparse.h"
20#include "colvardeps.h"
21
22#ifdef LEPTON
23#include "Lepton.h" // for runtime custom expressions
24#endif
25
50
51class colvar : public colvarparse, public colvardeps {
52
53public:
54
56 std::string name;
57
59 colvarvalue const & value() const;
60
62 colvarvalue const & actual_value() const;
63
65 colvarvalue const & run_ave() const;
66
68 cvm::real const & force_constant() const;
69
71 colvarvalue const & velocity() const;
72
81 colvarvalue const & total_force() const;
82
92
94 static std::vector<feature *> cv_features;
95
97 virtual const std::vector<feature *> &features() const
98 {
99 return cv_features;
100 }
101 virtual std::vector<feature *> &modify_features()
102 {
103 return cv_features;
104 }
105 static void delete_features() {
106 for (size_t i=0; i < cv_features.size(); i++) {
107 delete cv_features[i];
108 }
109 cv_features.clear();
110 }
111
115 void do_feature_side_effects(int id);
116
118 std::vector<colvarbias *> biases;
119
120protected:
121
122
123 /*
124 extended:
125 H = H_{0} + \sum_{i} 1/2*\lambda*(S_i(x(t))-s_i(t))^2 \\
126 + \sum_{i} 1/2*m_i*(ds_i(t)/dt)^2 \\
127 + \sum_{t'<t} W * exp(-1/2*\sum_{i} (s_i(t')-s_i(t))^2/(\delta{}s_i)^2) \\
128 + \sum_{w} (\sum_{i}c_{w,i}s_i(t) - D_w)^M/(\Sigma_w)^M
129
130 normal:
131 H = H_{0} + \sum_{t'<t} W * exp(-1/2*\sum_{i} (S_i(x(t'))-S_i(x(t)))^2/(\delta{}S_i)^2) \\
132 + \sum_{w} (\sum_{i}c_{w,i}S_i(t) - D_w)^M/(\Sigma_w)^M
133
134 output:
135 H = H_{0} (only output S(x), no forces)
136
137 Here:
138 S(x(t)) = x
139 s(t) = x_ext
140 DS = Ds = delta
141 */
142
143
146
147 // TODO: implement functionality to treat these
148 // /// Vector of individual values from CVCs
149 // colvarvalue x_cvc;
150
151 // /// Jacobian matrix of individual values from CVCs
152 // colvarvalue dx_cvc;
153
156
159
160 inline colvarvalue fdiff_velocity(colvarvalue const &xold,
161 colvarvalue const &xnew)
162 {
163 // using the gradient of the square distance to calculate the
164 // velocity (non-scalar variables automatically taken into
165 // account)
166 cvm::real dt = cvm::dt();
167 return ( ( (dt > 0.0) ? (1.0/dt) : 1.0 ) *
168 0.5 * dist2_lgrad(xnew, xold) );
169 }
170
173
174 // Options for extended_lagrangian
191
194
197
200
201public:
202
203
207
210
215
218
222
225
228
232
235
238
240 bool periodic_boundaries() const;
241
243 bool periodic_boundaries(colvarvalue const &lb, colvarvalue const &ub) const;
244
245
247 colvar();
248
250 int init(std::string const &conf);
251
253 int init_components(std::string const &conf);
254
256 int init_custom_function(std::string const &conf);
257
259 int init_grid_parameters(std::string const &conf);
260
262 int init_extended_Lagrangian(std::string const &conf);
263
265 int init_output_flags(std::string const &conf);
266
268 virtual int init_dependencies();
269
270private:
271
273 template<typename def_class_name> int init_components_type(std::string const & conf,
274 char const *def_desc,
275 char const *def_config_key);
276
279 int init_components_type_from_global_map(const std::string& conf,
280 const char* def_config_key);
281
282public:
283
285 void setup();
286
288 ~colvar();
289
290
292 int calc();
293
295 int end_of_step();
296
300 int calc_cvcs(int first = 0, size_t num_cvcs = 0);
301
303 int check_cvc_range(int first_cvc, size_t num_cvcs);
304
306 int calc_cvc_values(int first, size_t num_cvcs);
308 int calc_cvc_gradients(int first, size_t num_cvcs);
310 int calc_cvc_total_force(int first, size_t num_cvcs);
312 int calc_cvc_Jacobians(int first, size_t num_cvcs);
313
315 int collect_cvc_data();
316
318 int collect_cvc_values();
327
329 inline colvarvalue const applied_force() const
330 {
331 if (is_enabled(f_cv_extended_Lagrangian)) {
332 return fr;
333 }
334 return f;
335 }
336
338 void reset_bias_force();
339
341 void add_bias_force(colvarvalue const &force);
342
344 void add_bias_force_actual_value(colvarvalue const &force);
345
351
354
357 void communicate_forces();
358
360 int set_cvc_flags(std::vector<bool> const &flags);
361
363 int update_cvc_flags();
364
366 int update_cvc_config(std::vector<std::string> const &confs);
367
369 int cvc_param_exists(std::string const &param_name);
370
372 cvm::real get_cvc_param(std::string const &param_name);
373
375 void const *get_cvc_param_ptr(std::string const &param_name);
376
378 colvarvalue const *get_cvc_param_grad(std::string const &param_name);
379
381 int set_cvc_param(std::string const &param_name, void const *new_value);
382
383protected:
384
387
390
393
396
397public:
398
400 inline size_t num_dimensions() const
401 {
402 return value().size();
403 }
404
406 inline size_t num_cvcs() const
407 {
408 return cvcs.size();
409 }
410
413 inline size_t num_active_cvcs() const
414 {
415 return n_active_cvcs;
416 }
417
422 cvm::real dist2(colvarvalue const &x1,
423 colvarvalue const &x2) const;
424
430 colvarvalue const &x2) const;
431
437 colvarvalue const &x2) const;
438
443 void wrap(colvarvalue &x_unwrapped) const;
444
446 int parse_analysis(std::string const &conf);
447
449 int analyze();
450
452 std::istream & read_traj(std::istream &is);
453
455 std::ostream & write_traj(std::ostream &os);
457 std::ostream & write_traj_label(std::ostream &os);
458
460 std::istream & read_state(std::istream &is);
462 std::ostream & write_state(std::ostream &os);
463
465 int write_output_files();
466
467protected:
468
471
474
477
480 std::list< std::list<colvarvalue> > acf_x_history, acf_v_history;
483 std::list< std::list<colvarvalue> >::iterator acf_x_history_p, acf_v_history_p;
484
486 std::list< std::list<colvarvalue> > x_history;
489 std::list< std::list<colvarvalue> >::iterator x_history_p;
490
493 std::string acf_colvar_name;
505 std::vector<cvm::real> acf;
507 std::string acf_outfile;
508
520 };
521
524
526 void calc_vel_acf(std::list<colvarvalue> &v_history,
527 colvarvalue const &v);
528
531 void calc_coor_acf(std::list<colvarvalue> &x_history,
532 colvarvalue const &x);
533
536 void calc_p2coor_acf(std::list<colvarvalue> &x_history,
537 colvarvalue const &x);
538
540 int calc_acf();
542 int write_acf(std::ostream &os);
543
549 std::string runave_outfile;
554
556 int calc_runave();
557
562
563public:
564
565 // collective variable component base class
566 class cvc;
567
568 // list of available collective variable components
569
570 // scalar colvar components
571 class distance;
572 class distance_z;
573 class distance_xy;
574 class polar_theta;
575 class polar_phi;
576 class distance_inv;
577 class distance_pairs;
578 class dipole_magnitude;
579 class angle;
580 class dipole_angle;
581 class dihedral;
582 class coordnum;
583 class selfcoordnum;
584 class groupcoordnum;
585 class h_bond;
586 class rmsd;
587 class orientation_angle;
588 class orientation_proj;
589 class tilt;
590 class spin_angle;
591 class gyration;
592 class inertia;
593 class inertia_z;
594 class eigenvector;
595 class alpha_dihedrals;
596 class alpha_angles;
597 class dihedPC;
598 class alch_lambda;
599 class alch_Flambda;
600 class componentDisabled;
601 class CartesianBasedPath;
602 class gspath;
603 class gzpath;
604 class linearCombination;
605 class CVBasedPath;
606 class gspathCV;
607 class gzpathCV;
608 class aspathCV;
609 class azpathCV;
610 class euler_phi;
611 class euler_psi;
612 class euler_theta;
613 class neuralNetwork;
614 class customColvar;
615
616 // non-scalar components
617 class distance_vec;
618 class distance_dir;
619 class cartesian;
620 class orientation;
621
622 // components that do not handle any atoms directly
623 class map_total;
624
626 static const std::map<std::string, std::function<colvar::cvc *(const std::string &subcv_conf)>> &
628 {
629 return global_cvc_map;
630 }
631
633 static bool compare_cvc(const colvar::cvc* const i, const colvar::cvc* const j);
634
635protected:
636
638 std::vector<cvc *> cvcs;
639
641 std::vector<bool> cvc_flags;
642
645 void build_atom_list(void);
646
648 std::string scripted_function;
649
652 std::vector<const colvarvalue *> sorted_cvc_values;
653
654#ifdef LEPTON
656 std::vector<Lepton::CompiledExpression *> value_evaluators;
657
659 std::vector<Lepton::CompiledExpression *> gradient_evaluators;
660
662 std::vector<double *> value_eval_var_refs;
663 std::vector<double *> grad_eval_var_refs;
664
666 double dev_null;
667#endif
668
670 static std::map<std::string, std::function<colvar::cvc *(const std::string &subcv_conf)>>
672
674 std::vector<int> volmap_ids_;
675
676public:
677
679 std::vector<int> atom_ids;
680
684 std::vector<cvm::rvector> atomic_gradients;
685
687 virtual std::vector<std::vector<int> > get_atom_lists();
688
690 std::vector<int> const &get_volmap_ids();
691
692};
693
694
695inline cvm::real const & colvar::force_constant() const
696{
697 return ext_force_k;
698}
699
700
701inline colvarvalue const & colvar::value() const
702{
703 return x_reported;
704}
705
706
707inline colvarvalue const & colvar::actual_value() const
708{
709 return x;
710}
711
712
713inline colvarvalue const & colvar::run_ave() const
714{
715 return runave;
716}
717
718
719inline colvarvalue const & colvar::velocity() const
720{
721 return v_reported;
722}
723
724
725inline colvarvalue const & colvar::total_force() const
726{
727 return ft_reported;
728}
729
730
731inline void colvar::add_bias_force(colvarvalue const &force)
732{
734 std::string("applying a force to the variable \""+name+"\""));
735 if (cvm::debug()) {
736 cvm::log("Adding biasing force "+cvm::to_str(force)+" to colvar \""+name+"\".\n");
737 }
738 fb += force;
739}
740
741
743{
744 if (cvm::debug()) {
745 cvm::log("Adding biasing force "+cvm::to_str(force)+" to colvar \""+name+"\".\n");
746 }
747 fb_actual += force;
748}
749
750
752 fb.type(value());
753 fb.reset();
756}
757
758#endif
759
Definition: colvarcomp.h:1651
Definition: colvarcomp.h:1540
Definition: colvarcomp.h:1504
Definition: colvarcomp.h:1481
Colvar component: alpha helix content of a contiguous segment of 5 or more residues,...
Definition: colvarcomp.h:1129
Colvar component: angle between the centers of mass of three groups (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:748
Definition: colvarcomp.h:1712
Definition: colvarcomp.h:1726
Definition: colvarcomp.h:1461
Definition: colvarcomp.h:1525
Colvar component: coordination number between two groups (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:892
custom expression of colvars
Definition: colvarcomp.h:1626
Colvar component (base class for collective variables)
Definition: colvarcomp.h:73
Colvar component: dihedPC Projection of the config onto a dihedral principal component See e....
Definition: colvarcomp.h:1173
Colvar component: dihedral between the centers of mass of four groups (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:841
Colvar component: angle between the dipole of a molecule and an axis formed by two groups of atoms(co...
Definition: colvarcomp.h:795
Colvar component: dipole magnitude of a molecule.
Definition: colvarcomp.h:604
Colvar component: distance unit vector (direction) between centers of mass of two groups (colvarvalue...
Definition: colvarcomp.h:405
Colvar component: average distance between two groups of atoms, weighted as the sixth power,...
Definition: colvarcomp.h:556
Colvar component: N1xN2 vector of pairwise distances (colvarvalue::type_vector type,...
Definition: colvarcomp.h:584
Definition: colvarcomp.h:379
Colvar component: projection of the distance vector on a plane (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:471
Colvar component: projection of the distance vector along an axis(colvarvalue::type_scalar type,...
Definition: colvarcomp.h:430
Colvar component: distance between the centers of mass of two groups (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:347
Colvar component: projection of 3N coordinates onto an eigenvector(colvarvalue::type_scalar type,...
Definition: colvarcomp.h:709
Definition: colvarcomp.h:1352
Definition: colvarcomp.h:1374
Definition: colvarcomp.h:1396
Colvar component: coordination number between two groups (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:1015
Colvar component: alternative path collective variable using geometry, variable s Allow any combinati...
Definition: colvarcomp.h:1680
Colvar component: alternative path collective variable using geometry, variable s For more informatio...
Definition: colvarcomp.h:1566
Colvar component: Radius of gyration of an atom group (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:634
Definition: colvarcomp.h:1696
Colvar component: alternative path collective variable using geometry, variable z This should be merg...
Definition: colvarcomp.h:1587
Colvar component: hydrogen bond, defined as the product of a colvar::coordnum and 1/2*(1-cos((180-ang...
Definition: colvarcomp.h:1051
Colvar component: moment of inertia of an atom group around a user-defined axis (colvarvalue::type_sc...
Definition: colvarcomp.h:683
Colvar component: moment of inertia of an atom group (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:660
Current only linear combination of sub-CVCs is available.
Definition: colvarcomp.h:1606
Definition: colvarcomp.h:1764
Definition: colvarcomp.h:1745
Colvar component: angle of rotation with respect to a set of reference coordinates (colvarvalue::type...
Definition: colvarcomp.h:1247
Colvar component: cosine of the angle of rotation with respect to a set of reference coordinates (col...
Definition: colvarcomp.h:1271
Colvar component: orientation in space of an atom group, with respect to a set of reference coordinat...
Definition: colvarcomp.h:1205
Colvar component: polar coordinate phi of a group (colvarvalue::type_scalar type, range [-180:180])
Definition: colvarcomp.h:499
Colvar component: polar coordinate theta of a group (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:529
Colvar component: root mean square deviation (RMSD) of a group with respect to a set of reference coo...
Definition: colvarcomp.h:1422
Colvar component: self-coordination number within a group (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:975
Colvar component: angle of rotation around a predefined axis (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:1322
Colvar component: projection of the orientation vector onto a predefined axis (colvarvalue::type_scal...
Definition: colvarcomp.h:1295
A collective variable (main class); to be defined, it needs at least one object of a derived class of...
Definition: colvar.h:51
bool periodic_boundaries() const
Is the interval defined by the two boundaries periodic?
Definition: colvar.cpp:2157
void update_active_cvc_square_norm()
Update the sum of square coefficients for active cvcs.
Definition: colvar.cpp:2018
void communicate_forces()
Communicate forces (previously calculated in colvar::update()) to the external degrees of freedom.
Definition: colvar.cpp:1919
size_t num_cvcs() const
Number of CVC objects defined.
Definition: colvar.h:406
int end_of_step()
Carry out operations needed before next step is run.
Definition: colvar.cpp:1900
colvarvalue const * get_cvc_param_grad(std::string const &param_name)
Pointer to the derivative of the variable with respect to param_name.
Definition: colvar.cpp:2120
colvarvalue const & actual_value() const
Current actual value (not extended DOF)
Definition: colvar.h:707
colvarvalue const & run_ave() const
Current running average (if calculated as set by analysis flag)
Definition: colvar.h:713
colvarvalue const & velocity() const
Current velocity (previously set by calc() or by read_traj())
Definition: colvar.h:719
cvm::real active_cvc_square_norm
Sum of square coefficients for active cvcs.
Definition: colvar.h:389
colvarvalue x_restart
Value read from the most recent state file (if any)
Definition: colvar.h:473
size_t num_dimensions() const
Number of dimensions of the value of this colvar.
Definition: colvar.h:400
cvm::real runave_variance
Current value of the square deviation from the running average.
Definition: colvar.h:553
int calc_cvc_total_force(int first, size_t num_cvcs)
Same as colvar::calc_cvc_values but for total forces.
Definition: colvar.cpp:1573
colvarvalue f
Total applied force; fr (if extended_lagrangian is defined), fb (if biases are applied) and the walls...
Definition: colvar.h:214
std::list< std::list< colvarvalue > >::iterator x_history_p
Definition: colvar.h:489
std::vector< colvarbias * > biases
List of biases that depend on this colvar.
Definition: colvar.h:118
std::ostream & write_state(std::ostream &os)
Write the collective variable to a restart file.
Definition: colvar.cpp:2376
std::vector< const colvarvalue * > sorted_cvc_values
Definition: colvar.h:652
colvarvalue const & value() const
Current value (previously set by calc() or by read_traj())
Definition: colvar.h:701
int write_output_files()
Write output files (if defined, e.g. in analysis mode)
Definition: colvar.cpp:2518
int parse_analysis(std::string const &conf)
Read the analysis tasks.
Definition: colvar.cpp:1002
int init_custom_function(std::string const &conf)
Parse parameters for custom function with Lepton.
Definition: colvar.cpp:465
cvm::real ext_gamma
Friction coefficient for Langevin extended dynamics.
Definition: colvar.h:188
std::string acf_outfile
Name of the file to write the ACF.
Definition: colvar.h:507
cvm::real ext_mass
Mass of the restraint center.
Definition: colvar.h:184
std::ostream & write_traj(std::ostream &os)
Output formatted values to the trajectory file.
Definition: colvar.cpp:2464
cvm::real dist2(colvarvalue const &x1, colvarvalue const &x2) const
Use the internal metrics (as from colvar::cvc objects) to calculate square distances and gradients.
Definition: colvar.cpp:2168
bool after_restart
True if a state file was just read.
Definition: colvar.h:476
std::string name
Name.
Definition: colvar.h:56
colvarvalue v_ext
Velocity of the restraint center.
Definition: colvar.h:180
virtual std::vector< std::vector< int > > get_atom_lists()
Get vector of vectors of atom IDs for all atom groups.
Definition: colvar.cpp:1238
size_t acf_nframes
Number of frames for each ACF point.
Definition: colvar.h:501
int calc()
Calculate the colvar's value and related quantities.
Definition: colvar.cpp:1329
cvm::real period
Period, if this variable is periodic.
Definition: colvar.h:224
cvm::real kinetic_energy
If extended Lagrangian active: colvar kinetic energy.
Definition: colvar.h:559
void calc_vel_acf(std::list< colvarvalue > &v_history, colvarvalue const &v)
Velocity ACF, scalar product between v(0) and v(t)
Definition: colvar.cpp:2677
colvarvalue lower_boundary
Location of the lower boundary.
Definition: colvar.h:234
int update_cvc_config(std::vector< std::string > const &confs)
Modify the configuration of CVCs (currently, only base class data)
Definition: colvar.cpp:2053
static std::map< std::string, std::function< colvar::cvc *(const std::string &subcv_conf)> > global_cvc_map
A global mapping of cvc names to the cvc constructors.
Definition: colvar.h:671
int init_output_flags(std::string const &conf)
Init output flags.
Definition: colvar.cpp:707
int collect_cvc_values()
Collect the values of the CVCs.
Definition: colvar.cpp:1449
std::istream & read_traj(std::istream &is)
Read the value from a collective variable trajectory file.
Definition: colvar.cpp:2326
int check_cvc_range(int first_cvc, size_t num_cvcs)
Ensure that the selected range of CVCs is consistent.
Definition: colvar.cpp:1407
colvarvalue fb
Bias force; reset_bias_force() should be called before the biases are updated.
Definition: colvar.h:206
size_t n_active_cvcs
Number of CVC objects with an active flag.
Definition: colvar.h:386
std::vector< cvm::real > acf
ACF values.
Definition: colvar.h:505
int cvc_param_exists(std::string const &param_name)
Whether this named parameter exists (in the first and only component)
Definition: colvar.cpp:2088
std::list< std::list< colvarvalue > >::iterator acf_x_history_p
Definition: colvar.h:483
colvarvalue dist2_rgrad(colvarvalue const &x1, colvarvalue const &x2) const
Use the internal metrics (as from colvar::cvc objects) to calculate square distances and gradients.
Definition: colvar.cpp:2206
size_t runave_length
Length of running average series.
Definition: colvar.h:545
acf_type_e
Type of autocorrelation function (ACF)
Definition: colvar.h:510
@ acf_coor
Coordinate ACF, scalar product between x(0) and x(t)
Definition: colvar.h:516
@ acf_notset
Unset type.
Definition: colvar.h:512
@ acf_vel
Velocity ACF, scalar product between v(0) and v(t)
Definition: colvar.h:514
@ acf_p2coor
Coordinate ACF, second order Legendre polynomial between x(0) and x(t) (does not work with scalar num...
Definition: colvar.h:519
colvarvalue fj
Jacobian force, when Jacobian_force is enabled.
Definition: colvar.h:196
colvarvalue const applied_force() const
Get the current applied force.
Definition: colvar.h:329
void calc_p2coor_acf(std::list< colvarvalue > &x_history, colvarvalue const &x)
Coordinate ACF, second order Legendre polynomial between x(0) and x(t) (does not work with scalar num...
Definition: colvar.cpp:2722
int init_grid_parameters(std::string const &conf)
Init defaults for grid options.
Definition: colvar.cpp:481
int collect_cvc_Jacobians()
Same as colvar::collect_cvc_values but for Jacobian derivatives/forces.
Definition: colvar.cpp:1655
void wrap(colvarvalue &x_unwrapped) const
Use the internal metrics (as from colvar::cvc objects) to wrap a value into a standard interval.
Definition: colvar.cpp:2226
void do_feature_side_effects(int id)
Definition: colvar.cpp:945
int calc_runave()
Calculate the running average and its standard deviation.
Definition: colvar.cpp:2802
cvm::step_number prev_timestep
Absolute timestep number when this colvar was last updated.
Definition: colvar.h:395
colvarvalue prev_x_ext
Previous value of the restraint center;.
Definition: colvar.h:178
std::string scripted_function
Name of scripted function to be used.
Definition: colvar.h:648
void add_bias_force(colvarvalue const &force)
Add to the total force from biases.
Definition: colvar.h:731
virtual int init_dependencies()
Initialize dependency tree.
Definition: colvar.cpp:1075
std::ostream & write_traj_label(std::ostream &os)
Write a label to the trajectory file (comment line)
Definition: colvar.cpp:2413
int calc_colvar_properties()
Calculate the quantities associated to the colvar (but not to the CVCs)
Definition: colvar.cpp:1677
colvar()
Constructor.
Definition: colvar.cpp:29
int calc_cvc_Jacobians(int first, size_t num_cvcs)
Same as colvar::calc_cvc_values but for Jacobian derivatives/forces.
Definition: colvar.cpp:1634
static const std::map< std::string, std::function< colvar::cvc *(const std::string &subcv_conf)> > & get_global_cvc_map()
A global mapping of cvc names to the cvc constructors.
Definition: colvar.h:627
int collect_cvc_total_forces()
Same as colvar::collect_cvc_values but for total forces.
Definition: colvar.cpp:1602
bool expand_boundaries
Expand the boundaries of multiples of width, to keep the value always within range.
Definition: colvar.h:231
void calc_coor_acf(std::list< colvarvalue > &x_history, colvarvalue const &x)
Coordinate ACF, scalar product between x(0) and x(t) (does not work with scalar numbers)
Definition: colvar.cpp:2702
static bool compare_cvc(const colvar::cvc *const i, const colvar::cvc *const j)
function for sorting cvcs by their names
Definition: colvar.cpp:49
std::list< std::list< colvarvalue > > x_history
Time series of values and velocities used in running averages.
Definition: colvar.h:486
colvarvalue v_reported
Cached reported velocity.
Definition: colvar.h:172
colvarvalue fr
Applied force on extended DOF, for output (unscaled if using MTS)
Definition: colvar.h:193
colvarvalue upper_boundary
Location of the upper boundary.
Definition: colvar.h:237
std::list< std::list< colvarvalue > > acf_x_history
Definition: colvar.h:480
bool acf_normalize
Normalize the ACF to a maximum value of 1?
Definition: colvar.h:503
colvarvalue v_fdiff
Finite-difference velocity.
Definition: colvar.h:158
std::vector< cvc * > cvcs
Array of colvar::cvc objects.
Definition: colvar.h:638
colvarvalue prev_v_ext
Previous velocity of the restraint center.
Definition: colvar.h:182
int init_extended_Lagrangian(std::string const &conf)
Init extended Lagrangian parameters.
Definition: colvar.cpp:629
static std::vector< feature * > cv_features
Implementation of the feature list for colvar.
Definition: colvar.h:94
cvm::real ext_sigma
Amplitude of Gaussian white noise for Langevin extended dynamics.
Definition: colvar.h:190
std::vector< cvm::rvector > atomic_gradients
Array of atomic gradients collected from all cvcs with appropriate components, rotations etc....
Definition: colvar.h:684
cvm::real width
Typical fluctuation amplitude for this collective variable (e.g. local width of a free energy basin)
Definition: colvar.h:91
int init_components_type(std::string const &conf, char const *def_desc, char const *def_config_key)
Parse the CVC configuration for all components of a certain type.
Definition: colvar.cpp:743
colvarvalue runave
Current value of the running average.
Definition: colvar.h:551
void const * get_cvc_param_ptr(std::string const &param_name)
Get a pointer to the named parameter (from the first and only component)
Definition: colvar.cpp:2109
size_t acf_offset
After how many steps the ACF starts.
Definition: colvar.h:497
colvarvalue ft_reported
Cached reported total force.
Definition: colvar.h:199
acf_type_e acf_type
Type of autocorrelation function (ACF)
Definition: colvar.h:523
void build_atom_list(void)
Initialize the sorted list of atom IDs for atoms involved in all cvcs (called when enabling f_cv_coll...
Definition: colvar.cpp:962
std::istream & read_state(std::istream &is)
Read the collective variable from a restart file.
Definition: colvar.cpp:2245
cvm::real get_cvc_param(std::string const &param_name)
Get the value of the named parameter (from the first and only component)
Definition: colvar.cpp:2098
cvm::real wrap_center
Center of wrapping, if this variable is periodic.
Definition: colvar.h:227
void add_bias_force_actual_value(colvarvalue const &force)
Apply a force to the actual value (only meaningful with extended Lagrangian)
Definition: colvar.h:742
colvarvalue f_old
Applied force at the previous step (to be subtracted from total force if needed)
Definition: colvar.h:217
int update_cvc_flags()
Updates the flags in the CVC objects, and their number.
Definition: colvar.cpp:2029
cvm::real const & force_constant() const
Force constant of the spring.
Definition: colvar.h:695
cvm::real ext_force_k
Restraint force constant.
Definition: colvar.h:186
colvarvalue dist2_lgrad(colvarvalue const &x1, colvarvalue const &x2) const
Use the internal metrics (as from colvar::cvc objects) to calculate square distances and gradients.
Definition: colvar.cpp:2187
~colvar()
Destructor.
Definition: colvar.cpp:1264
int calc_acf()
Calculate the auto-correlation function (ACF)
Definition: colvar.cpp:2582
void reset_bias_force()
Set the total biasing force to zero.
Definition: colvar.h:751
colvarvalue ft
Total force, as derived from the atomic trajectory; should equal the system force plus f.
Definition: colvar.h:221
std::string acf_colvar_name
Collective variable with which the correlation is calculated (default: itself)
Definition: colvar.h:493
colvarvalue const & total_force() const
Current total force (previously obtained from calc() or read_traj()). Note: this is calculated using ...
Definition: colvar.h:725
int init_components(std::string const &conf)
Parse the CVC configuration and allocate their data.
Definition: colvar.cpp:830
int collect_cvc_gradients()
Same as colvar::collect_cvc_values but for gradients.
Definition: colvar.cpp:1556
colvarvalue fb_actual
Bias force to the actual value (only useful with extended Lagrangian)
Definition: colvar.h:209
int collect_cvc_data()
Collect quantities from CVCs and update aggregated data for the colvar.
Definition: colvar.cpp:1378
size_t runave_stride
Timesteps to skip between two values in the running average series.
Definition: colvar.h:547
int calc_cvcs(int first=0, size_t num_cvcs=0)
Calculate a subset of the colvar components (CVCs) currently active (default: all active CVCs) Note: ...
Definition: colvar.cpp:1344
colvarvalue x_reported
Cached reported value (x may be manipulated)
Definition: colvar.h:155
colvarvalue x_old
Previous value (to calculate velocities during analysis)
Definition: colvar.h:470
int init(std::string const &conf)
Main init function.
Definition: colvar.cpp:55
cvm::real update_forces_energy()
Collect all forces on this colvar, integrate internal equations of motion of internal degrees of free...
Definition: colvar.cpp:1743
std::vector< int > const & get_volmap_ids()
Volmap numeric IDs, one for each CVC (-1 if not available)
Definition: colvar.cpp:1249
cvm::real potential_energy
If extended Lagrangian active: colvar harmonic potential.
Definition: colvar.h:561
colvarvalue x
Value of the colvar.
Definition: colvar.h:145
size_t acf_stride
How many timesteps separate two ACF values.
Definition: colvar.h:499
size_t num_active_cvcs() const
number of CVC objects with an active flag (as set by update_cvc_flags)
Definition: colvar.h:413
void update_extended_Lagrangian()
Integrate equations of motion of extended Lagrangian coordinate if needed.
Definition: colvar.cpp:1789
int calc_cvc_values(int first, size_t num_cvcs)
Calculate the values of the given subset of CVCs.
Definition: colvar.cpp:1418
int init_components_type_from_global_map(const std::string &conf, const char *def_config_key)
Definition: colvar.cpp:755
std::vector< bool > cvc_flags
Flags to enable or disable cvcs at next colvar evaluation.
Definition: colvar.h:641
std::vector< int > volmap_ids_
Volmap numeric IDs, one for each CVC (-1 if not available)
Definition: colvar.h:674
size_t acf_length
Length of autocorrelation function (ACF)
Definition: colvar.h:495
std::string runave_outfile
Name of the file to write the running average.
Definition: colvar.h:549
int set_cvc_param(std::string const &param_name, void const *new_value)
Set the named parameter in the first and only component to the given value.
Definition: colvar.cpp:2131
int write_acf(std::ostream &os)
Save the ACF to a file.
Definition: colvar.cpp:2744
int calc_cvc_gradients(int first, size_t num_cvcs)
Same as colvar::calc_cvc_values but for gradients.
Definition: colvar.cpp:1521
virtual const std::vector< feature * > & features() const
Implementation of the feature list accessor for colvar.
Definition: colvar.h:97
colvarvalue x_ext
Restraint center.
Definition: colvar.h:176
int set_cvc_flags(std::vector< bool > const &flags)
Enables and disables individual CVCs based on the given array.
Definition: colvar.cpp:2005
int analyze()
Perform analysis tasks.
Definition: colvar.cpp:2548
std::vector< int > atom_ids
Sorted array of (zero-based) IDs for all atoms involved.
Definition: colvar.h:679
void setup()
Get ready for a run and re-initialize internal data if needed.
Definition: colvar.cpp:1224
Parent class for a member object of a bias, cv or cvc etc. containing features and their dependencies...
Definition: colvardeps.h:34
void check_enabled(int f, std::string const &reason) const
Check that a feature is enabled, raising COLVARS_BUG_ERROR if not.
Definition: colvardeps.h:430
@ f_cv_extended_Lagrangian
The variable has a harmonic restraint around a moving center with fictitious mass; bias forces will b...
Definition: colvardeps.h:289
@ f_cv_gradient
Gradients are calculated and temporarily stored, so that external forces can be applied.
Definition: colvardeps.h:266
double real
Defining an abstract real number allows to switch precision.
Definition: colvarmodule.h:100
static real dt()
Time step of MD integrator (fs)
Definition: colvarmodule.cpp:2068
static void log(std::string const &message, int min_log_level=10)
Definition: colvarmodule.cpp:1722
static bool debug()
Whether debug output should be enabled (compile-time option)
Definition: colvarmodule.h:333
static std::string to_str(char const *s)
Convert to string for output purposes.
Definition: colvarmodule.cpp:2136
long long step_number
Use a 64-bit integer to store the step number.
Definition: colvarmodule.h:97
Base class containing parsing functions; all objects which need to parse input inherit from this.
Definition: colvarparse.h:26
Value of a collective variable: this is a metatype which can be set at runtime. By default it is set ...
Definition: colvarvalue.h:41
Type type() const
Get the current type.
Definition: colvarvalue.h:154
void reset()
Set to the null value for the data type currently defined.
Definition: colvarvalue.cpp:202
size_t size() const
Number of dimensions of this variable.
Definition: colvarvalue.h:354
Collective variables main module.
Parsing functions for collective variables.