Collective Variables Module - Developer Documentation
Loading...
Searching...
No Matches
colvarbias_meta.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 COLVARBIAS_META_H
11#define COLVARBIAS_META_H
12
13#include <iosfwd>
14#include <list>
15#include <memory>
16#include <vector>
17
18#include "colvarbias.h"
19
22
23
24
27 : public virtual colvarbias,
28 public virtual colvarbias_ti
29{
30
31public:
32
39 };
40
43
44 colvarbias_meta(char const *key);
45 virtual ~colvarbias_meta();
46
47 virtual int init(std::string const &conf);
48 virtual int init_replicas_params(std::string const &conf);
49 virtual int init_well_tempered_params(std::string const &conf);
50 virtual int init_ebmeta_params(std::string const &conf);
51
52 virtual int clear_state_data();
53
54 virtual int update();
55 virtual int update_grid_params();
56 virtual int update_bias();
57 virtual int update_grid_data();
58 virtual int replica_share();
59 virtual size_t replica_share_freq() const;
60
61 virtual int calc_energy(std::vector<colvarvalue> const *values);
62 virtual int calc_forces(std::vector<colvarvalue> const *values);
63
64 virtual std::string const get_state_params() const;
65 virtual int set_state_params(std::string const &state_conf);
66
67 virtual std::ostream &write_state_data(std::ostream &os);
69 virtual std::istream &read_state_data(std::istream &is);
71
72private:
73
74 template <typename IST, typename GT>
75 IST &read_grid_data_template_(IST &is, std::string const &key, GT *grid, GT *backup_grid);
76
77 template <typename IST> IST &read_state_data_template_(IST &is);
78
79 template <typename OST> OST &write_state_data_template_(OST &os);
80
81public:
82
85
86 virtual int setup_output();
87 virtual int write_output_files();
88 virtual void write_pmf();
89 virtual int write_state_to_replicas();
90
91 class hill;
92 typedef std::list<hill>::iterator hill_iter;
93
94protected:
95
101
103 std::vector<cvm::real> colvar_sigmas;
104
107
110
112 std::string const hills_traj_file_name() const;
113
116 std::list<hill> hills;
117
121
124 std::list<hill> hills_off_grid;
125
128
130 void recount_hills_off_grid(hill_iter h_first, hill_iter h_last);
131
132 template <typename OST> OST &write_hill_template_(OST &os, colvarbias_meta::hill const &h);
133
135 std::ostream &write_hill(std::ostream &os, hill const &h);
136
139
140 template <typename IST> IST &read_hill_template_(IST &is);
141
143 std::istream & read_hill(std::istream &is);
144
147
151 std::list<hill>::const_iterator add_hill(hill const &h);
152
155 std::list<hill>::const_iterator delete_hill(hill_iter &h);
156
159 virtual void calc_hills(hill_iter h_first,
160 hill_iter h_last,
161 cvm::real &energy,
162 std::vector<colvarvalue> const *values);
163
167 virtual void calc_hills_force(size_t const &i,
168 hill_iter h_first,
169 hill_iter h_last,
170 std::vector<colvarvalue> &forces,
171 std::vector<colvarvalue> const *values);
172
173
176
180
183
186
189
192
195
198
202
206
209
212
214 bool ebmeta;
215
217 std::unique_ptr<colvar_grid_scalar> target_dist;
218
221
222
227
229 std::shared_ptr<colvar_grid_scalar> hills_energy;
230
232 std::shared_ptr<colvar_grid_gradient> hills_energy_gradients;
233
235 void project_hills(hill_iter h_first, hill_iter h_last, colvar_grid_scalar *ge,
236 colvar_grid_gradient *gf, bool print_progress = false);
237
238
239 // Multiple Replicas variables and functions
240
242 std::string replica_id;
243
245 std::string replica_file_name;
246
248 virtual int update_replicas_registry();
249
251 virtual int read_replica_files();
252
254 virtual int write_replica_state_file();
255
257 virtual int reopen_replica_buffer_file();
258
264 std::vector<colvarbias_meta *> replicas;
265
268
273 std::string replicas_registry;
275 std::string replica_list_file;
276
282
286
291
294
296 std::ostringstream hills_traj_os_buf;
297};
298
299
300
301
304
305protected:
306
309
312
315
318
320 std::vector<colvarvalue> centers;
321
323 std::vector<cvm::real> sigmas;
324
326 std::string replica;
327
328public:
329
330 friend class colvarbias_meta;
331
339 std::vector<colvarvalue> const &cv_values,
340 std::vector<cvm::real> const &cv_sigmas,
341 std::string const &replica = "");
342
344 hill(colvarbias_meta::hill const &h);
345
347 ~hill();
348
351
354 {
355 return W * sW * hill_value;
356 }
357
359 inline cvm::real energy(cvm::real const &new_weight)
360 {
361 return new_weight * sW * hill_value;
362 }
363
365 inline cvm::real const &value()
366 {
367 return hill_value;
368 }
369
371 inline void value(cvm::real const &new_value)
372 {
373 hill_value = new_value;
374 }
375
378 {
379 return W * sW;
380 }
381
383 inline void scale(cvm::real const &new_scale_fac)
384 {
385 sW = new_scale_fac;
386 }
387
389 inline std::vector<colvarvalue> & center()
390 {
391 return centers;
392 }
393
395 inline colvarvalue & center(size_t const &i)
396 {
397 return centers[i];
398 }
399
401 inline friend bool operator < (hill const &h1, hill const &h2)
402 {
403 if (h1.it < h2.it) return true;
404 else return false;
405 }
406
408 inline friend bool operator <= (hill const &h1, hill const &h2)
409 {
410 if (h1.it <= h2.it) return true;
411 else return false;
412 }
413
415 inline friend bool operator > (hill const &h1, hill const &h2)
416 {
417 if (h1.it > h2.it) return true;
418 else return false;
419 }
420
422 inline friend bool operator >= (hill const &h1, hill const &h2)
423 {
424 if (h1.it >= h2.it) return true;
425 else return false;
426 }
427
429 inline friend bool operator == (hill const &h1, hill const &h2)
430 {
431 if ( (h1.it >= h2.it) && (h1.replica == h2.replica) ) return true;
432 else return false;
433 }
434
436 std::string output_traj();
437
438};
439
440
441#endif
Class for accumulating the gradient of a scalar function on a grid.
Definition: colvargrid.h:1598
Class for accumulating a scalar function on a grid.
Definition: colvargrid.h:1283
A hill for the metadynamics bias.
Definition: colvarbias_meta.h:303
colvarvalue & center(size_t const &i)
Get the i-th component of the center.
Definition: colvarbias_meta.h:395
friend bool operator>=(hill const &h1, hill const &h2)
Comparison operator.
Definition: colvarbias_meta.h:422
cvm::real energy(cvm::real const &new_weight)
Get the energy using another hill weight.
Definition: colvarbias_meta.h:359
cvm::real hill_value
Value of the hill function (ranges between 0 and 1)
Definition: colvarbias_meta.h:311
~hill()
Destructor.
Definition: colvarbias_meta.cpp:2111
cvm::step_number it
Time step at which this hill was added.
Definition: colvarbias_meta.h:308
void value(cvm::real const &new_value)
Set the hill value as specified.
Definition: colvarbias_meta.h:371
hill & operator=(colvarbias_meta::hill const &h)
Assignment operator.
Definition: colvarbias_meta.cpp:2097
cvm::real const & value()
Get the current hill value.
Definition: colvarbias_meta.h:365
friend bool operator<=(hill const &h1, hill const &h2)
Comparison operator.
Definition: colvarbias_meta.h:408
friend bool operator<(hill const &h1, hill const &h2)
Comparison operator.
Definition: colvarbias_meta.h:401
std::string replica
Identity of the replica who added this hill.
Definition: colvarbias_meta.h:326
std::vector< colvarvalue > centers
Centers of the hill in the collective variable space.
Definition: colvarbias_meta.h:320
cvm::real W
Maximum height in energy of the hill.
Definition: colvarbias_meta.h:317
cvm::real sW
Scale factor, which could be modified at runtime (default: 1)
Definition: colvarbias_meta.h:314
friend bool operator==(hill const &h1, hill const &h2)
Comparison operator.
Definition: colvarbias_meta.h:429
void scale(cvm::real const &new_scale_fac)
Scale the weight with this factor (by default 1.0 is used)
Definition: colvarbias_meta.h:383
cvm::real weight()
Get the weight.
Definition: colvarbias_meta.h:377
std::vector< colvarvalue > & center()
Get the center of the hill.
Definition: colvarbias_meta.h:389
std::string output_traj()
Represent the hill ina string suitable for a trajectory file.
Definition: colvarbias_meta.cpp:2024
friend bool operator>(hill const &h1, hill const &h2)
Comparison operator.
Definition: colvarbias_meta.h:415
std::vector< cvm::real > sigmas
Half-widths of the hill in the collective variable space.
Definition: colvarbias_meta.h:323
cvm::real energy()
Get the energy.
Definition: colvarbias_meta.h:353
Metadynamics bias (implementation of colvarbias)
Definition: colvarbias_meta.h:29
bool dump_replica_fes
Dump the free energy surface (.pmf file) every restartFrequency using only the hills from this replic...
Definition: colvarbias_meta.h:201
bool safely_read_restart
Try to read the restart information by allocating new grids before replacing the current ones (used e...
Definition: colvarbias_meta.h:226
std::string replica_id
Identifier for this replica.
Definition: colvarbias_meta.h:242
virtual int reopen_replica_buffer_file()
Call this after write_replica_state_file()
Definition: colvarbias_meta.cpp:2006
std::list< hill > hills_off_grid
List of hills used on this bias that are on the boundary edges; these are updated regardless of wheth...
Definition: colvarbias_meta.h:124
size_t grids_freq
How often the hills should be projected onto the grids.
Definition: colvarbias_meta.h:188
virtual int update()
Definition: colvarbias_meta.cpp:402
cvm::real hill_weight
Height of new hills.
Definition: colvarbias_meta.h:175
bool use_grids
Bin the hills on grids of energy and forces, and use them to force the colvars (as opposed to derivin...
Definition: colvarbias_meta.h:179
size_t new_hill_freq
Number of simulation steps between two hills.
Definition: colvarbias_meta.h:106
virtual void calc_hills_force(size_t const &i, hill_iter h_first, hill_iter h_last, std::vector< colvarvalue > &forces, std::vector< colvarvalue > const *values)
Calculate the forces acting on the i-th colvar, incrementing colvar_forces[i]; must be called after c...
Definition: colvarbias_meta.cpp:788
std::string replicas_registry_file
Definition: colvarbias_meta.h:271
virtual int calc_forces(std::vector< colvarvalue > const *values)
Definition: colvarbias_meta.cpp:688
std::unique_ptr< colvar_grid_scalar > target_dist
Target distribution for EBmeta.
Definition: colvarbias_meta.h:217
bool ebmeta
Ensemble-biased metadynamics (EBmeta) flag.
Definition: colvarbias_meta.h:214
virtual int write_output_files()
Write any output files that this bias may have (e.g. PMF files)
Definition: colvarbias_meta.cpp:1877
std::string replica_list_file
List of files written by this replica.
Definition: colvarbias_meta.h:275
std::ostream & write_hill(std::ostream &os, hill const &h)
Write a hill to a formatted stream.
Definition: colvarbias_meta.cpp:1548
std::shared_ptr< colvar_grid_scalar > hills_energy
Hill energy, cached on a grid.
Definition: colvarbias_meta.h:229
virtual std::ostream & write_state_data(std::ostream &os)
Write all mutable data not already written by get_state_params() to a formatted stream.
Definition: colvarbias_meta.cpp:1850
cvm::step_number ebmeta_equil_steps
Number of equilibration steps for EBmeta.
Definition: colvarbias_meta.h:220
void project_hills(hill_iter h_first, hill_iter h_last, colvar_grid_scalar *ge, colvar_grid_gradient *gf, bool print_progress=false)
Project the selected hills onto grids.
Definition: colvarbias_meta.cpp:863
std::list< hill > hills
List of hills used on this bias (total); if a grid is employed, these don't need to be updated at eve...
Definition: colvarbias_meta.h:116
virtual int init(std::string const &conf)
Parse config string and (re)initialize.
Definition: colvarbias_meta.cpp:50
bool keep_hills
Keep hills in the restart file (e.g. to accurately rebin later)
Definition: colvarbias_meta.h:191
std::string replica_hills_file
Definition: colvarbias_meta.h:290
virtual std::string const get_state_params() const
Write the values of specific mutable properties to a string.
Definition: colvarbias_meta.cpp:1803
virtual int clear_state_data()
Delete only the allocatable data (save memory)
Definition: colvarbias_meta.cpp:324
virtual std::istream & read_state_data(std::istream &is)
Read all mutable data not already set by set_state_params() from a formatted stream.
Definition: colvarbias_meta.cpp:1423
cvm::real bias_temperature
Biasing temperature in well-tempered metadynamics.
Definition: colvarbias_meta.h:211
bool rebin_grids
Rebin the hills upon restarting.
Definition: colvarbias_meta.h:182
virtual int read_replica_files()
Read new data from replicas' files.
Definition: colvarbias_meta.cpp:1121
std::list< hill >::const_iterator add_hill(hill const &h)
Add a new hill; if a .hills trajectory is written, write it there; if there is more than one replica,...
Definition: colvarbias_meta.cpp:337
std::string replica_state_file
Definition: colvarbias_meta.h:279
std::istream & read_hill(std::istream &is)
Read a new hill from a formatted stream.
Definition: colvarbias_meta.cpp:1684
bool expand_grids
Should the grids be expanded if necessary?
Definition: colvarbias_meta.h:185
bool well_tempered
Whether to use well-tempered metadynamics.
Definition: colvarbias_meta.h:208
bool restart_keep_hills
value of keepHills saved in the most recent restart file
Definition: colvarbias_meta.h:194
std::ostringstream hills_traj_os_buf
Cache of the hills trajectory.
Definition: colvarbias_meta.h:296
bool dump_fes_save
Dump the free energy surface files at different time steps, appending the step number to each file.
Definition: colvarbias_meta.h:205
virtual int setup_output()
(Re)initialize the output files (does not write them yet)
Definition: colvarbias_meta.cpp:1696
size_t update_status
Definition: colvarbias_meta.h:285
Communication
Communication between different replicas.
Definition: colvarbias_meta.h:34
@ multiple_replicas
Hills added concurrently by several replicas.
Definition: colvarbias_meta.h:38
@ single_replica
One replica (default)
Definition: colvarbias_meta.h:36
size_t replica_update_freq
Frequency at which data the "mirror" biases are updated.
Definition: colvarbias_meta.h:267
Communication comm
Communication between different replicas.
Definition: colvarbias_meta.h:42
virtual int write_replica_state_file()
Write full state information to be read by other replicas.
Definition: colvarbias_meta.cpp:1978
bool b_hills_traj
Write the hill logfile.
Definition: colvarbias_meta.h:109
virtual size_t replica_share_freq() const
Report the frequency at which this bias needs to communicate with replicas.
Definition: colvarbias_meta.cpp:977
virtual void calc_hills(hill_iter h_first, hill_iter h_last, cvm::real &energy, std::vector< colvarvalue > const *values)
Calculate the values of the hills, incrementing bias_energy.
Definition: colvarbias_meta.cpp:758
std::shared_ptr< colvar_grid_gradient > hills_energy_gradients
Hill forces, cached on a grid.
Definition: colvarbias_meta.h:232
cvm::real hill_width
Definition: colvarbias_meta.h:100
bool replica_state_file_in_sync
Whether a mirror bias has read the latest version of its state file.
Definition: colvarbias_meta.h:281
void rebin_grids_after_restart()
Function called by read_state_data() to execute rebinning (if requested)
Definition: colvarbias_meta.cpp:1435
std::streampos replica_hills_file_pos
Position within replica_hills_file (when reading it)
Definition: colvarbias_meta.h:293
std::list< hill >::const_iterator delete_hill(hill_iter &h)
Remove a previously saved hill (returns an iterator for the next hill in the list)
Definition: colvarbias_meta.cpp:370
virtual int write_state_to_replicas()
If this bias is communicating with other replicas through files, send it to them.
Definition: colvarbias_meta.cpp:1862
void recount_hills_off_grid(hill_iter h_first, hill_iter h_last)
Regenerate the hills_off_grid list.
Definition: colvarbias_meta.cpp:941
virtual int calc_energy(std::vector< colvarvalue > const *values)
Definition: colvarbias_meta.cpp:627
virtual int set_state_params(std::string const &state_conf)
Read the values of specific mutable properties from a string.
Definition: colvarbias_meta.cpp:1250
hill_iter new_hills_begin
Iterator to the first of the "newest" hills (when using grids, those who haven't been mapped yet)
Definition: colvarbias_meta.h:120
std::string replica_file_name
File containing the paths to the output files from this replica.
Definition: colvarbias_meta.h:245
std::string const hills_traj_file_name() const
Name of the hill logfile.
Definition: colvarbias_meta.cpp:1792
virtual int update_replicas_registry()
Read the existing replicas on registry.
Definition: colvarbias_meta.cpp:983
std::string replicas_registry
List of replicas (and their output list files)
Definition: colvarbias_meta.h:273
bool dump_fes
Dump the free energy surface (.pmf file) every restartFrequency.
Definition: colvarbias_meta.h:197
hill_iter new_hills_off_grid_begin
Same as new_hills_begin, but for the off-grid ones.
Definition: colvarbias_meta.h:127
std::vector< colvarbias_meta * > replicas
Additional, "mirror" metadynamics biases, to collect info from the other replicas.
Definition: colvarbias_meta.h:264
std::vector< cvm::real > colvar_sigmas
The sigma parameters of the Gaussian hills.
Definition: colvarbias_meta.h:103
Base class for unconstrained thermodynamic-integration FE estimator.
Definition: colvarbias.h:333
Collective variable bias, base class.
Definition: colvarbias.h:23
double real
Defining an abstract real number allows to switch precision.
Definition: colvarmodule.h:100
long long step_number
Use a 64-bit integer to store the step number.
Definition: colvarmodule.h:97
Value of a collective variable: this is a metatype which can be set at runtime. By default it is set ...
Definition: colvarvalue.h:43
Definition: colvars_memstream.h:30