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
15#if (__cplusplus >= 201103L)
16#include <map>
17#include <functional>
18#endif
19
20#include "colvarmodule.h"
21#include "colvarvalue.h"
22#include "colvarparse.h"
23#include "colvardeps.h"
24
25#ifdef LEPTON
26#include "Lepton.h" // for runtime custom expressions
27#endif
28
53
54class colvar : public colvarparse, public colvardeps {
55
56public:
57
59 std::string name;
60
62 colvarvalue const & value() const;
63
65 colvarvalue const & actual_value() const;
66
68 colvarvalue const & run_ave() const;
69
71 cvm::real const & force_constant() const;
72
74 colvarvalue const & velocity() const;
75
84 colvarvalue const & total_force() const;
85
95
97 static std::vector<feature *> cv_features;
98
100 virtual const std::vector<feature *> &features() const
101 {
102 return cv_features;
103 }
104 virtual std::vector<feature *> &modify_features()
105 {
106 return cv_features;
107 }
108 static void delete_features() {
109 for (size_t i=0; i < cv_features.size(); i++) {
110 delete cv_features[i];
111 }
112 cv_features.clear();
113 }
114
118 void do_feature_side_effects(int id);
119
121 std::vector<colvarbias *> biases;
122
123protected:
124
125
126 /*
127 extended:
128 H = H_{0} + \sum_{i} 1/2*\lambda*(S_i(x(t))-s_i(t))^2 \\
129 + \sum_{i} 1/2*m_i*(ds_i(t)/dt)^2 \\
130 + \sum_{t'<t} W * exp(-1/2*\sum_{i} (s_i(t')-s_i(t))^2/(\delta{}s_i)^2) \\
131 + \sum_{w} (\sum_{i}c_{w,i}s_i(t) - D_w)^M/(\Sigma_w)^M
132
133 normal:
134 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) \\
135 + \sum_{w} (\sum_{i}c_{w,i}S_i(t) - D_w)^M/(\Sigma_w)^M
136
137 output:
138 H = H_{0} (only output S(x), no forces)
139
140 Here:
141 S(x(t)) = x
142 s(t) = x_ext
143 DS = Ds = delta
144 */
145
146
149
150 // TODO: implement functionality to treat these
151 // /// Vector of individual values from CVCs
152 // colvarvalue x_cvc;
153
154 // /// Jacobian matrix of individual values from CVCs
155 // colvarvalue dx_cvc;
156
159
162
163 inline colvarvalue fdiff_velocity(colvarvalue const &xold,
164 colvarvalue const &xnew)
165 {
166 // using the gradient of the square distance to calculate the
167 // velocity (non-scalar variables automatically taken into
168 // account)
169 cvm::real dt = cvm::dt();
170 return ( ( (dt > 0.0) ? (1.0/dt) : 1.0 ) *
171 0.5 * dist2_lgrad(xnew, xold) );
172 }
173
176
177 // Options for extended_lagrangian
194
197
200
203
204public:
205
206
210
213
218
221
225
228
231
235
238
241
243 bool periodic_boundaries() const;
244
246 bool periodic_boundaries(colvarvalue const &lb, colvarvalue const &ub) const;
247
248
250 colvar();
251
253 int init(std::string const &conf);
254
256 int init_components(std::string const &conf);
257
259 int init_custom_function(std::string const &conf);
260
262 int init_grid_parameters(std::string const &conf);
263
265 int init_extended_Lagrangian(std::string const &conf);
266
268 int init_output_flags(std::string const &conf);
269
271 virtual int init_dependencies();
272
273private:
275 template<typename def_class_name> int init_components_type(std::string const & conf,
276 char const *def_desc,
277 char const *def_config_key);
278#if (__cplusplus >= 201103L)
282 int init_components_type_from_global_map(const std::string& conf,
283 const char* def_config_key);
284#endif
285
286public:
287
289 void setup();
290
292 ~colvar();
293
294
296 int calc();
297
299 int end_of_step();
300
304 int calc_cvcs(int first = 0, size_t num_cvcs = 0);
305
307 int check_cvc_range(int first_cvc, size_t num_cvcs);
308
310 int calc_cvc_values(int first, size_t num_cvcs);
312 int calc_cvc_gradients(int first, size_t num_cvcs);
314 int calc_cvc_total_force(int first, size_t num_cvcs);
316 int calc_cvc_Jacobians(int first, size_t num_cvcs);
317
319 int collect_cvc_data();
320
322 int collect_cvc_values();
331
333 inline colvarvalue const applied_force() const
334 {
335 if (is_enabled(f_cv_extended_Lagrangian)) {
336 return fr;
337 }
338 return f;
339 }
340
342 void reset_bias_force();
343
345 void add_bias_force(colvarvalue const &force);
346
348 void add_bias_force_actual_value(colvarvalue const &force);
349
355
358
361 void communicate_forces();
362
364 int set_cvc_flags(std::vector<bool> const &flags);
365
367 int update_cvc_flags();
368
370 int update_cvc_config(std::vector<std::string> const &confs);
371
373 int cvc_param_exists(std::string const &param_name);
374
376 cvm::real get_cvc_param(std::string const &param_name);
377
379 void const *get_cvc_param_ptr(std::string const &param_name);
380
382 colvarvalue const *get_cvc_param_grad(std::string const &param_name);
383
385 int set_cvc_param(std::string const &param_name, void const *new_value);
386
387protected:
388
391
394
397
400
401public:
402
404 inline size_t num_dimensions() const
405 {
406 return value().size();
407 }
408
410 inline size_t num_cvcs() const
411 {
412 return cvcs.size();
413 }
414
417 inline size_t num_active_cvcs() const
418 {
419 return n_active_cvcs;
420 }
421
426 cvm::real dist2(colvarvalue const &x1,
427 colvarvalue const &x2) const;
428
434 colvarvalue const &x2) const;
435
441 colvarvalue const &x2) const;
442
447 void wrap(colvarvalue &x_unwrapped) const;
448
450 int parse_analysis(std::string const &conf);
451
453 int analyze();
454
456 std::istream & read_traj(std::istream &is);
457
459 std::ostream & write_traj(std::ostream &os);
461 std::ostream & write_traj_label(std::ostream &os);
462
464 std::istream & read_state(std::istream &is);
466 std::ostream & write_state(std::ostream &os);
467
469 int write_output_files();
470
471protected:
472
475
478
481
484 std::list< std::list<colvarvalue> > acf_x_history, acf_v_history;
487 std::list< std::list<colvarvalue> >::iterator acf_x_history_p, acf_v_history_p;
488
490 std::list< std::list<colvarvalue> > x_history;
493 std::list< std::list<colvarvalue> >::iterator x_history_p;
494
497 std::string acf_colvar_name;
509 std::vector<cvm::real> acf;
511 std::string acf_outfile;
512
524 };
525
528
530 void calc_vel_acf(std::list<colvarvalue> &v_history,
531 colvarvalue const &v);
532
535 void calc_coor_acf(std::list<colvarvalue> &x_history,
536 colvarvalue const &x);
537
540 void calc_p2coor_acf(std::list<colvarvalue> &x_history,
541 colvarvalue const &x);
542
544 int calc_acf();
546 int write_acf(std::ostream &os);
547
553 std::string runave_outfile;
558
560 int calc_runave();
561
566
567public:
568
569 // collective variable component base class
570 class cvc;
571
572 // list of available collective variable components
573
574 // scalar colvar components
575 class distance;
576 class distance_z;
577 class distance_xy;
578 class polar_theta;
579 class polar_phi;
580 class distance_inv;
581 class distance_pairs;
582 class dipole_magnitude;
583 class angle;
584 class dipole_angle;
585 class dihedral;
586 class coordnum;
587 class selfcoordnum;
588 class groupcoordnum;
589 class h_bond;
590 class rmsd;
591 class orientation_angle;
592 class orientation_proj;
593 class tilt;
594 class spin_angle;
595 class gyration;
596 class inertia;
597 class inertia_z;
598 class eigenvector;
599 class alpha_dihedrals;
600 class alpha_angles;
601 class dihedPC;
602 class alch_lambda;
603 class alch_Flambda;
604 class componentDisabled;
605 class CartesianBasedPath;
606 class gspath;
607 class gzpath;
608 class linearCombination;
609 class CVBasedPath;
610 class gspathCV;
611 class gzpathCV;
612 class aspathCV;
613 class azpathCV;
614 class euler_phi;
615 class euler_psi;
616 class euler_theta;
617 class neuralNetwork;
618 class customColvar;
619
620 // non-scalar components
621 class distance_vec;
622 class distance_dir;
623 class cartesian;
624 class orientation;
625
626 // components that do not handle any atoms directly
627 class map_total;
628
630#if (__cplusplus >= 201103L)
632 static const std::map<std::string, std::function<colvar::cvc* (const std::string& subcv_conf)>>& get_global_cvc_map() {
633 return global_cvc_map;
634 }
635#endif
636
638 static bool compare_cvc(const colvar::cvc* const i, const colvar::cvc* const j);
639
640protected:
641
643 std::vector<cvc *> cvcs;
644
646 std::vector<bool> cvc_flags;
647
650 void build_atom_list(void);
651
653 std::string scripted_function;
654
657 std::vector<const colvarvalue *> sorted_cvc_values;
658
659#ifdef LEPTON
661 std::vector<Lepton::CompiledExpression *> value_evaluators;
662
664 std::vector<Lepton::CompiledExpression *> gradient_evaluators;
665
667 std::vector<double *> value_eval_var_refs;
668 std::vector<double *> grad_eval_var_refs;
669
671 double dev_null;
672#endif
673
674#if (__cplusplus >= 201103L)
676 static std::map<std::string, std::function<colvar::cvc* (const std::string& subcv_conf)>> global_cvc_map;
677#endif
678
680 std::vector<int> volmap_ids_;
681
682public:
683
685 std::vector<int> atom_ids;
686
690 std::vector<cvm::rvector> atomic_gradients;
691
693 virtual std::vector<std::vector<int> > get_atom_lists();
694
696 std::vector<int> const &get_volmap_ids();
697
698};
699
700
701inline cvm::real const & colvar::force_constant() const
702{
703 return ext_force_k;
704}
705
706
707inline colvarvalue const & colvar::value() const
708{
709 return x_reported;
710}
711
712
713inline colvarvalue const & colvar::actual_value() const
714{
715 return x;
716}
717
718
719inline colvarvalue const & colvar::run_ave() const
720{
721 return runave;
722}
723
724
725inline colvarvalue const & colvar::velocity() const
726{
727 return v_reported;
728}
729
730
731inline colvarvalue const & colvar::total_force() const
732{
733 return ft_reported;
734}
735
736
737inline void colvar::add_bias_force(colvarvalue const &force)
738{
740 std::string("applying a force to the variable \""+name+"\""));
741 if (cvm::debug()) {
742 cvm::log("Adding biasing force "+cvm::to_str(force)+" to colvar \""+name+"\".\n");
743 }
744 fb += force;
745}
746
747
749{
750 if (cvm::debug()) {
751 cvm::log("Adding biasing force "+cvm::to_str(force)+" to colvar \""+name+"\".\n");
752 }
753 fb_actual += force;
754}
755
756
758 fb.type(value());
759 fb.reset();
762}
763
764#endif
765
Definition: colvarcomp.h:1782
Definition: colvarcomp.h:1775
Definition: colvarcomp.h:1509
Definition: colvarcomp.h:1486
Colvar component: alpha helix content of a contiguous segment of 5 or more residues,...
Definition: colvarcomp.h:1134
Colvar component: angle between the centers of mass of three groups (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:753
Definition: colvarcomp.h:1817
Definition: colvarcomp.h:1824
Definition: colvarcomp.h:1466
Definition: colvarcomp.h:1530
Colvar component: coordination number between two groups (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:897
Colvar component (base class for collective variables)
Definition: colvarcomp.h:78
Colvar component: dihedPC Projection of the config onto a dihedral principal component See e....
Definition: colvarcomp.h:1178
Colvar component: dihedral between the centers of mass of four groups (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:846
Colvar component: angle between the dipole of a molecule and an axis formed by two groups of atoms(co...
Definition: colvarcomp.h:800
Colvar component: dipole magnitude of a molecule.
Definition: colvarcomp.h:609
Colvar component: distance unit vector (direction) between centers of mass of two groups (colvarvalue...
Definition: colvarcomp.h:410
Colvar component: average distance between two groups of atoms, weighted as the sixth power,...
Definition: colvarcomp.h:561
Colvar component: N1xN2 vector of pairwise distances (colvarvalue::type_vector type,...
Definition: colvarcomp.h:589
Definition: colvarcomp.h:384
Colvar component: projection of the distance vector on a plane (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:476
Colvar component: projection of the distance vector along an axis(colvarvalue::type_scalar type,...
Definition: colvarcomp.h:435
Colvar component: distance between the centers of mass of two groups (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:352
Colvar component: projection of 3N coordinates onto an eigenvector(colvarvalue::type_scalar type,...
Definition: colvarcomp.h:714
Definition: colvarcomp.h:1357
Definition: colvarcomp.h:1379
Definition: colvarcomp.h:1401
Colvar component: coordination number between two groups (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:1020
Definition: colvarcomp.h:1803
Definition: colvarcomp.h:1789
Colvar component: Radius of gyration of an atom group (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:639
Definition: colvarcomp.h:1810
Definition: colvarcomp.h:1796
Colvar component: hydrogen bond, defined as the product of a colvar::coordnum and 1/2*(1-cos((180-ang...
Definition: colvarcomp.h:1056
Colvar component: moment of inertia of an atom group around a user-defined axis (colvarvalue::type_sc...
Definition: colvarcomp.h:688
Colvar component: moment of inertia of an atom group (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:665
Definition: colvarcomp.h:1768
Definition: colvarcomp.h:1843
Definition: colvarcomp.h:1831
Colvar component: angle of rotation with respect to a set of reference coordinates (colvarvalue::type...
Definition: colvarcomp.h:1252
Colvar component: cosine of the angle of rotation with respect to a set of reference coordinates (col...
Definition: colvarcomp.h:1276
Colvar component: orientation in space of an atom group, with respect to a set of reference coordinat...
Definition: colvarcomp.h:1210
Colvar component: polar coordinate phi of a group (colvarvalue::type_scalar type, range [-180:180])
Definition: colvarcomp.h:504
Colvar component: polar coordinate theta of a group (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:534
Colvar component: root mean square deviation (RMSD) of a group with respect to a set of reference coo...
Definition: colvarcomp.h:1427
Colvar component: self-coordination number within a group (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:980
Colvar component: angle of rotation around a predefined axis (colvarvalue::type_scalar type,...
Definition: colvarcomp.h:1327
Colvar component: projection of the orientation vector onto a predefined axis (colvarvalue::type_scal...
Definition: colvarcomp.h:1300
A collective variable (main class); to be defined, it needs at least one object of a derived class of...
Definition: colvar.h:54
bool periodic_boundaries() const
Is the interval defined by the two boundaries periodic?
Definition: colvar.cpp:2168
void update_active_cvc_square_norm()
Update the sum of square coefficients for active cvcs.
Definition: colvar.cpp:2029
void communicate_forces()
Communicate forces (previously calculated in colvar::update()) to the external degrees of freedom.
Definition: colvar.cpp:1930
size_t num_cvcs() const
Number of CVC objects defined.
Definition: colvar.h:410
int end_of_step()
Carry out operations needed before next step is run.
Definition: colvar.cpp:1911
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:2131
colvarvalue const & actual_value() const
Current actual value (not extended DOF)
Definition: colvar.h:713
colvarvalue const & run_ave() const
Current running average (if calculated as set by analysis flag)
Definition: colvar.h:719
colvarvalue const & velocity() const
Current velocity (previously set by calc() or by read_traj())
Definition: colvar.h:725
cvm::real active_cvc_square_norm
Sum of square coefficients for active cvcs.
Definition: colvar.h:393
colvarvalue x_restart
Value read from the most recent state file (if any)
Definition: colvar.h:477
size_t num_dimensions() const
Number of dimensions of the value of this colvar.
Definition: colvar.h:404
cvm::real runave_variance
Current value of the square deviation from the running average.
Definition: colvar.h:557
int calc_cvc_total_force(int first, size_t num_cvcs)
Same as colvar::calc_cvc_values but for total forces.
Definition: colvar.cpp:1584
colvarvalue f
Total applied force; fr (if extended_lagrangian is defined), fb (if biases are applied) and the walls...
Definition: colvar.h:217
std::list< std::list< colvarvalue > >::iterator x_history_p
Definition: colvar.h:493
std::vector< colvarbias * > biases
List of biases that depend on this colvar.
Definition: colvar.h:121
std::ostream & write_state(std::ostream &os)
Write the collective variable to a restart file.
Definition: colvar.cpp:2387
std::vector< const colvarvalue * > sorted_cvc_values
Definition: colvar.h:657
colvarvalue const & value() const
Current value (previously set by calc() or by read_traj())
Definition: colvar.h:707
int write_output_files()
Write output files (if defined, e.g. in analysis mode)
Definition: colvar.cpp:2529
int parse_analysis(std::string const &conf)
Read the analysis tasks.
Definition: colvar.cpp:1013
int init_custom_function(std::string const &conf)
Parse parameters for custom function with Lepton.
Definition: colvar.cpp:463
cvm::real ext_gamma
Friction coefficient for Langevin extended dynamics.
Definition: colvar.h:191
std::string acf_outfile
Name of the file to write the ACF.
Definition: colvar.h:511
cvm::real ext_mass
Mass of the restraint center.
Definition: colvar.h:187
std::ostream & write_traj(std::ostream &os)
Output formatted values to the trajectory file.
Definition: colvar.cpp:2475
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:2179
bool after_restart
True if a state file was just read.
Definition: colvar.h:480
std::string name
Name.
Definition: colvar.h:59
colvarvalue v_ext
Velocity of the restraint center.
Definition: colvar.h:183
virtual std::vector< std::vector< int > > get_atom_lists()
Get vector of vectors of atom IDs for all atom groups.
Definition: colvar.cpp:1249
size_t acf_nframes
Number of frames for each ACF point.
Definition: colvar.h:505
int calc()
Calculate the colvar's value and related quantities.
Definition: colvar.cpp:1340
cvm::real period
Period, if this variable is periodic.
Definition: colvar.h:227
cvm::real kinetic_energy
If extended Lagrangian active: colvar kinetic energy.
Definition: colvar.h:563
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:2688
colvarvalue lower_boundary
Location of the lower boundary.
Definition: colvar.h:237
int update_cvc_config(std::vector< std::string > const &confs)
Modify the configuration of CVCs (currently, only base class data)
Definition: colvar.cpp:2064
int init_output_flags(std::string const &conf)
Init output flags.
Definition: colvar.cpp:712
int collect_cvc_values()
Collect the values of the CVCs.
Definition: colvar.cpp:1460
std::istream & read_traj(std::istream &is)
Read the value from a collective variable trajectory file.
Definition: colvar.cpp:2337
int check_cvc_range(int first_cvc, size_t num_cvcs)
Ensure that the selected range of CVCs is consistent.
Definition: colvar.cpp:1418
colvarvalue fb
Bias force; reset_bias_force() should be called before the biases are updated.
Definition: colvar.h:209
size_t n_active_cvcs
Number of CVC objects with an active flag.
Definition: colvar.h:390
std::vector< cvm::real > acf
ACF values.
Definition: colvar.h:509
int cvc_param_exists(std::string const &param_name)
Whether this named parameter exists (in the first and only component)
Definition: colvar.cpp:2099
std::list< std::list< colvarvalue > >::iterator acf_x_history_p
Definition: colvar.h:487
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:2217
size_t runave_length
Length of running average series.
Definition: colvar.h:549
acf_type_e
Type of autocorrelation function (ACF)
Definition: colvar.h:514
@ acf_coor
Coordinate ACF, scalar product between x(0) and x(t)
Definition: colvar.h:520
@ acf_notset
Unset type.
Definition: colvar.h:516
@ acf_vel
Velocity ACF, scalar product between v(0) and v(t)
Definition: colvar.h:518
@ acf_p2coor
Coordinate ACF, second order Legendre polynomial between x(0) and x(t) (does not work with scalar num...
Definition: colvar.h:523
colvarvalue fj
Jacobian force, when Jacobian_force is enabled.
Definition: colvar.h:199
colvarvalue const applied_force() const
Get the current applied force.
Definition: colvar.h:333
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:2733
int init_grid_parameters(std::string const &conf)
Init defaults for grid options.
Definition: colvar.cpp:486
int collect_cvc_Jacobians()
Same as colvar::collect_cvc_values but for Jacobian derivatives/forces.
Definition: colvar.cpp:1666
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:2237
void do_feature_side_effects(int id)
Definition: colvar.cpp:956
int calc_runave()
Calculate the running average and its standard deviation.
Definition: colvar.cpp:2813
cvm::step_number prev_timestep
Absolute timestep number when this colvar was last updated.
Definition: colvar.h:399
colvarvalue prev_x_ext
Previous value of the restraint center;.
Definition: colvar.h:181
std::string scripted_function
Name of scripted function to be used.
Definition: colvar.h:653
void add_bias_force(colvarvalue const &force)
Add to the total force from biases.
Definition: colvar.h:737
virtual int init_dependencies()
Initialize dependency tree.
Definition: colvar.cpp:1086
std::ostream & write_traj_label(std::ostream &os)
Write a label to the trajectory file (comment line)
Definition: colvar.cpp:2424
int calc_colvar_properties()
Calculate the quantities associated to the colvar (but not to the CVCs)
Definition: colvar.cpp:1688
colvar()
Constructor.
Definition: colvar.cpp:27
int calc_cvc_Jacobians(int first, size_t num_cvcs)
Same as colvar::calc_cvc_values but for Jacobian derivatives/forces.
Definition: colvar.cpp:1645
int collect_cvc_total_forces()
Same as colvar::collect_cvc_values but for total forces.
Definition: colvar.cpp:1613
bool expand_boundaries
Expand the boundaries of multiples of width, to keep the value always within range.
Definition: colvar.h:234
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:2713
static bool compare_cvc(const colvar::cvc *const i, const colvar::cvc *const j)
getter of the global cvc map
Definition: colvar.cpp:47
std::list< std::list< colvarvalue > > x_history
Time series of values and velocities used in running averages.
Definition: colvar.h:490
colvarvalue v_reported
Cached reported velocity.
Definition: colvar.h:175
colvarvalue fr
Applied force on extended DOF, for output (unscaled if using MTS)
Definition: colvar.h:196
colvarvalue upper_boundary
Location of the upper boundary.
Definition: colvar.h:240
std::list< std::list< colvarvalue > > acf_x_history
Definition: colvar.h:484
bool acf_normalize
Normalize the ACF to a maximum value of 1?
Definition: colvar.h:507
colvarvalue v_fdiff
Finite-difference velocity.
Definition: colvar.h:161
std::vector< cvc * > cvcs
Array of colvar::cvc objects.
Definition: colvar.h:643
colvarvalue prev_v_ext
Previous velocity of the restraint center.
Definition: colvar.h:185
int init_extended_Lagrangian(std::string const &conf)
Init extended Lagrangian parameters.
Definition: colvar.cpp:634
static std::vector< feature * > cv_features
Implementation of the feature list for colvar.
Definition: colvar.h:97
cvm::real ext_sigma
Amplitude of Gaussian white noise for Langevin extended dynamics.
Definition: colvar.h:193
std::vector< cvm::rvector > atomic_gradients
Array of atomic gradients collected from all cvcs with appropriate components, rotations etc....
Definition: colvar.h:690
cvm::real width
Typical fluctuation amplitude for this collective variable (e.g. local width of a free energy basin)
Definition: colvar.h:94
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:760
colvarvalue runave
Current value of the running average.
Definition: colvar.h:555
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:2120
size_t acf_offset
After how many steps the ACF starts.
Definition: colvar.h:501
colvarvalue ft_reported
Cached reported total force.
Definition: colvar.h:202
acf_type_e acf_type
Type of autocorrelation function (ACF)
Definition: colvar.h:527
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:973
std::istream & read_state(std::istream &is)
Read the collective variable from a restart file.
Definition: colvar.cpp:2256
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:2109
cvm::real wrap_center
Center of wrapping, if this variable is periodic.
Definition: colvar.h:230
void add_bias_force_actual_value(colvarvalue const &force)
Apply a force to the actual value (only meaningful with extended Lagrangian)
Definition: colvar.h:748
colvarvalue f_old
Applied force at the previous step (to be subtracted from total force if needed)
Definition: colvar.h:220
int update_cvc_flags()
Updates the flags in the CVC objects, and their number.
Definition: colvar.cpp:2040
cvm::real const & force_constant() const
Force constant of the spring.
Definition: colvar.h:701
cvm::real ext_force_k
Restraint force constant.
Definition: colvar.h:189
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:2198
~colvar()
Destructor.
Definition: colvar.cpp:1275
int calc_acf()
Calculate the auto-correlation function (ACF)
Definition: colvar.cpp:2593
void reset_bias_force()
Set the total biasing force to zero.
Definition: colvar.h:757
colvarvalue ft
Total force, as derived from the atomic trajectory; should equal the system force plus f.
Definition: colvar.h:224
std::string acf_colvar_name
Collective variable with which the correlation is calculated (default: itself)
Definition: colvar.h:497
colvarvalue const & total_force() const
Current total force (previously obtained from calc() or read_traj()). Note: this is calculated using ...
Definition: colvar.h:731
int init_components(std::string const &conf)
Parse the CVC configuration and allocate their data.
Definition: colvar.cpp:841
int collect_cvc_gradients()
Same as colvar::collect_cvc_values but for gradients.
Definition: colvar.cpp:1567
colvarvalue fb_actual
Bias force to the actual value (only useful with extended Lagrangian)
Definition: colvar.h:212
int collect_cvc_data()
Collect quantities from CVCs and update aggregated data for the colvar.
Definition: colvar.cpp:1389
size_t runave_stride
Timesteps to skip between two values in the running average series.
Definition: colvar.h:551
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:1355
colvarvalue x_reported
Cached reported value (x may be manipulated)
Definition: colvar.h:158
colvarvalue x_old
Previous value (to calculate velocities during analysis)
Definition: colvar.h:474
int init(std::string const &conf)
Main init function.
Definition: colvar.cpp:53
cvm::real update_forces_energy()
Collect all forces on this colvar, integrate internal equations of motion of internal degrees of free...
Definition: colvar.cpp:1754
std::vector< int > const & get_volmap_ids()
Volmap numeric IDs, one for each CVC (-1 if not available)
Definition: colvar.cpp:1260
cvm::real potential_energy
If extended Lagrangian active: colvar harmonic potential.
Definition: colvar.h:565
colvarvalue x
Value of the colvar.
Definition: colvar.h:148
size_t acf_stride
How many timesteps separate two ACF values.
Definition: colvar.h:503
size_t num_active_cvcs() const
number of CVC objects with an active flag (as set by update_cvc_flags)
Definition: colvar.h:417
void update_extended_Lagrangian()
Integrate equations of motion of extended Lagrangian coordinate if needed.
Definition: colvar.cpp:1800
int calc_cvc_values(int first, size_t num_cvcs)
Calculate the values of the given subset of CVCs.
Definition: colvar.cpp:1429
std::vector< bool > cvc_flags
Flags to enable or disable cvcs at next colvar evaluation.
Definition: colvar.h:646
std::vector< int > volmap_ids_
Volmap numeric IDs, one for each CVC (-1 if not available)
Definition: colvar.h:680
size_t acf_length
Length of autocorrelation function (ACF)
Definition: colvar.h:499
std::string runave_outfile
Name of the file to write the running average.
Definition: colvar.h:553
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:2142
int write_acf(std::ostream &os)
Save the ACF to a file.
Definition: colvar.cpp:2755
int calc_cvc_gradients(int first, size_t num_cvcs)
Same as colvar::calc_cvc_values but for gradients.
Definition: colvar.cpp:1532
virtual const std::vector< feature * > & features() const
Implementation of the feature list accessor for colvar.
Definition: colvar.h:100
colvarvalue x_ext
Restraint center.
Definition: colvar.h:179
int set_cvc_flags(std::vector< bool > const &flags)
Enables and disables individual CVCs based on the given array.
Definition: colvar.cpp:2016
int analyze()
Perform analysis tasks.
Definition: colvar.cpp:2559
std::vector< int > atom_ids
Sorted array of (zero-based) IDs for all atoms involved.
Definition: colvar.h:685
void setup()
Get ready for a run and re-initialize internal data if needed.
Definition: colvar.cpp:1235
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:1734
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.