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);
46
47 int init(std::string const &conf) override;
48 int init_replicas_params(std::string const &conf);
49 int init_well_tempered_params(std::string const &conf);
50 int init_ebmeta_params(std::string const &conf);
51
52 int clear_state_data() override;
53
54 int update() override;
55 int update_grid_params();
56 int update_bias();
57 int update_grid_data();
58
59 int replica_share() override;
60 size_t replica_share_freq() const override;
61
62 int calc_energy(std::vector<colvarvalue> const *values) override;
63 int calc_forces(std::vector<colvarvalue> const *values) override;
64
65 std::string const get_state_params() const override;
66 int set_state_params(std::string const &state_conf) override;
67
68 std::ostream &write_state_data(std::ostream &os) override;
70 std::istream &read_state_data(std::istream &is) override;
72
73private:
74
75 template <typename IST, typename GT>
76 IST &read_grid_data_template_(IST &is, std::string const &key, GT *grid, GT *backup_grid);
77
78 template <typename IST> IST &read_state_data_template_(IST &is);
79
80 template <typename OST> OST &write_state_data_template_(OST &os);
81
82public:
83
86
87 int setup_output() override;
88 int write_output_files() override;
89 int write_state_to_replicas() override;
90 void write_pmf();
91
92 class hill;
93 typedef std::list<hill>::iterator hill_iter;
94
95protected:
96
102
104 std::vector<cvm::real> colvar_sigmas;
105
107 size_t new_hill_freq = 0;
108
111
113 std::string const hills_traj_file_name() const;
114
117 std::list<hill> hills;
118
122
125 std::list<hill> hills_off_grid;
126
129
131 void recount_hills_off_grid(hill_iter h_first, hill_iter h_last);
132
133 template <typename OST> OST &write_hill_template_(OST &os, colvarbias_meta::hill const &h);
134
136 std::ostream &write_hill(std::ostream &os, hill const &h);
137
140
141 template <typename IST> IST &read_hill_template_(IST &is);
142
144 std::istream & read_hill(std::istream &is);
145
148
152 std::list<hill>::const_iterator add_hill(hill const &h);
153
156 std::list<hill>::const_iterator delete_hill(hill_iter &h);
157
160 virtual void calc_hills(hill_iter h_first,
161 hill_iter h_last,
162 cvm::real &energy,
163 std::vector<colvarvalue> const *values);
164
168 virtual void calc_hills_force(size_t const &i,
169 hill_iter h_first,
170 hill_iter h_last,
171 std::vector<colvarvalue> &forces,
172 std::vector<colvarvalue> const *values);
173
174
177
181
184
187
189 size_t grids_freq = 0;
190
193
196
199
203
207
210
213
215 bool ebmeta;
216
218 std::unique_ptr<colvar_grid_scalar> target_dist;
219
222
223
228
230 std::shared_ptr<colvar_grid_scalar> hills_energy;
231
233 std::shared_ptr<colvar_grid_gradient> hills_energy_gradients;
234
236 void project_hills(hill_iter h_first, hill_iter h_last, colvar_grid_scalar *ge,
237 colvar_grid_gradient *gf, bool print_progress = false);
238
239
240 // Multiple Replicas variables and functions
241
243 std::string replica_id;
244
246 std::string replica_file_name;
247
249 virtual int update_replicas_registry();
250
252 virtual int read_replica_files();
253
255 virtual int write_replica_state_file();
256
258 virtual int reopen_replica_buffer_file();
259
265 std::vector<colvarbias_meta *> replicas;
266
269
274 std::string replicas_registry;
276 std::string replica_list_file;
277
283
287
292
295
297 std::ostringstream hills_traj_os_buf;
298};
299
300
301
302
305
306protected:
307
310
313
316
319
321 std::vector<colvarvalue> centers;
322
324 std::vector<cvm::real> sigmas;
325
327 std::string replica;
328
329public:
330
331 friend class colvarbias_meta;
332
340 std::vector<colvarvalue> const &cv_values,
341 std::vector<cvm::real> const &cv_sigmas,
342 std::string const &replica = "");
343
345 hill(colvarbias_meta::hill const &h);
346
348 ~hill();
349
352
355 {
356 return W * sW * hill_value;
357 }
358
360 inline cvm::real energy(cvm::real const &new_weight)
361 {
362 return new_weight * sW * hill_value;
363 }
364
366 inline cvm::real const &value()
367 {
368 return hill_value;
369 }
370
372 inline void value(cvm::real const &new_value)
373 {
374 hill_value = new_value;
375 }
376
379 {
380 return W * sW;
381 }
382
384 inline void scale(cvm::real const &new_scale_fac)
385 {
386 sW = new_scale_fac;
387 }
388
390 inline std::vector<colvarvalue> & center()
391 {
392 return centers;
393 }
394
396 inline colvarvalue & center(size_t const &i)
397 {
398 return centers[i];
399 }
400
402 inline friend bool operator < (hill const &h1, hill const &h2)
403 {
404 if (h1.it < h2.it) return true;
405 else return false;
406 }
407
409 inline friend bool operator <= (hill const &h1, hill const &h2)
410 {
411 if (h1.it <= h2.it) return true;
412 else return false;
413 }
414
416 inline friend bool operator > (hill const &h1, hill const &h2)
417 {
418 if (h1.it > h2.it) return true;
419 else return false;
420 }
421
423 inline friend bool operator >= (hill const &h1, hill const &h2)
424 {
425 if (h1.it >= h2.it) return true;
426 else return false;
427 }
428
430 inline friend bool operator == (hill const &h1, hill const &h2)
431 {
432 if ( (h1.it >= h2.it) && (h1.replica == h2.replica) ) return true;
433 else return false;
434 }
435
437 std::string output_traj();
438
439};
440
441
442#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:304
colvarvalue & center(size_t const &i)
Get the i-th component of the center.
Definition: colvarbias_meta.h:396
friend bool operator>=(hill const &h1, hill const &h2)
Comparison operator.
Definition: colvarbias_meta.h:423
cvm::real energy(cvm::real const &new_weight)
Get the energy using another hill weight.
Definition: colvarbias_meta.h:360
cvm::real hill_value
Value of the hill function (ranges between 0 and 1)
Definition: colvarbias_meta.h:312
~hill()
Destructor.
Definition: colvarbias_meta.cpp:2158
cvm::step_number it
Time step at which this hill was added.
Definition: colvarbias_meta.h:309
void value(cvm::real const &new_value)
Set the hill value as specified.
Definition: colvarbias_meta.h:372
hill & operator=(colvarbias_meta::hill const &h)
Assignment operator.
Definition: colvarbias_meta.cpp:2144
cvm::real const & value()
Get the current hill value.
Definition: colvarbias_meta.h:366
friend bool operator<=(hill const &h1, hill const &h2)
Comparison operator.
Definition: colvarbias_meta.h:409
friend bool operator<(hill const &h1, hill const &h2)
Comparison operator.
Definition: colvarbias_meta.h:402
std::string replica
Identity of the replica who added this hill.
Definition: colvarbias_meta.h:327
std::vector< colvarvalue > centers
Centers of the hill in the collective variable space.
Definition: colvarbias_meta.h:321
cvm::real W
Maximum height in energy of the hill.
Definition: colvarbias_meta.h:318
cvm::real sW
Scale factor, which could be modified at runtime (default: 1)
Definition: colvarbias_meta.h:315
friend bool operator==(hill const &h1, hill const &h2)
Comparison operator.
Definition: colvarbias_meta.h:430
void scale(cvm::real const &new_scale_fac)
Scale the weight with this factor (by default 1.0 is used)
Definition: colvarbias_meta.h:384
cvm::real weight()
Get the weight.
Definition: colvarbias_meta.h:378
std::vector< colvarvalue > & center()
Get the center of the hill.
Definition: colvarbias_meta.h:390
std::string output_traj()
Represent the hill ina string suitable for a trajectory file.
Definition: colvarbias_meta.cpp:2071
friend bool operator>(hill const &h1, hill const &h2)
Comparison operator.
Definition: colvarbias_meta.h:416
std::vector< cvm::real > sigmas
Half-widths of the hill in the collective variable space.
Definition: colvarbias_meta.h:324
cvm::real energy()
Get the energy.
Definition: colvarbias_meta.h:354
Metadynamics bias (implementation of colvarbias)
Definition: colvarbias_meta.h:29
std::string const get_state_params() const override
Write the values of specific mutable properties to a string.
Definition: colvarbias_meta.cpp:1850
bool dump_replica_fes
Dump the free energy surface (.pmf file) every restartFrequency using only the hills from this replic...
Definition: colvarbias_meta.h:202
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:227
std::string replica_id
Identifier for this replica.
Definition: colvarbias_meta.h:243
int write_output_files() override
Write any output files that this bias may have (e.g. PMF files)
Definition: colvarbias_meta.cpp:1924
virtual int reopen_replica_buffer_file()
Call this after write_replica_state_file()
Definition: colvarbias_meta.cpp:2053
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:125
int clear_state_data() override
Delete only the allocatable data (save memory)
Definition: colvarbias_meta.cpp:346
int set_state_params(std::string const &state_conf) override
Read the values of specific mutable properties from a string.
Definition: colvarbias_meta.cpp:1296
size_t grids_freq
How often the hills should be projected onto the grids.
Definition: colvarbias_meta.h:189
cvm::real hill_weight
Height of new hills.
Definition: colvarbias_meta.h:176
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:180
size_t new_hill_freq
Number of simulation steps between two hills.
Definition: colvarbias_meta.h:107
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:830
std::string replicas_registry_file
Definition: colvarbias_meta.h:272
int calc_energy(std::vector< colvarvalue > const *values) override
Definition: colvarbias_meta.cpp:669
std::unique_ptr< colvar_grid_scalar > target_dist
Target distribution for EBmeta.
Definition: colvarbias_meta.h:218
bool ebmeta
Ensemble-biased metadynamics (EBmeta) flag.
Definition: colvarbias_meta.h:215
std::string replica_list_file
List of files written by this replica.
Definition: colvarbias_meta.h:276
std::ostream & write_hill(std::ostream &os, hill const &h)
Write a hill to a formatted stream.
Definition: colvarbias_meta.cpp:1594
std::shared_ptr< colvar_grid_scalar > hills_energy
Hill energy, cached on a grid.
Definition: colvarbias_meta.h:230
cvm::step_number ebmeta_equil_steps
Number of equilibration steps for EBmeta.
Definition: colvarbias_meta.h:221
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:905
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:117
bool keep_hills
Keep hills in the restart file (e.g. to accurately rebin later)
Definition: colvarbias_meta.h:192
std::string replica_hills_file
Definition: colvarbias_meta.h:291
int calc_forces(std::vector< colvarvalue > const *values) override
Definition: colvarbias_meta.cpp:730
cvm::real bias_temperature
Biasing temperature in well-tempered metadynamics.
Definition: colvarbias_meta.h:212
bool rebin_grids
Rebin the hills upon restarting.
Definition: colvarbias_meta.h:183
virtual int read_replica_files()
Read new data from replicas' files.
Definition: colvarbias_meta.cpp:1163
int write_state_to_replicas() override
If this bias is communicating with other replicas through files, send it to them.
Definition: colvarbias_meta.cpp:1909
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:359
std::string replica_state_file
Definition: colvarbias_meta.h:280
std::istream & read_hill(std::istream &is)
Read a new hill from a formatted stream.
Definition: colvarbias_meta.cpp:1730
bool expand_grids
Should the grids be expanded if necessary?
Definition: colvarbias_meta.h:186
bool well_tempered
Whether to use well-tempered metadynamics.
Definition: colvarbias_meta.h:209
int update() override
Definition: colvarbias_meta.cpp:424
bool restart_keep_hills
value of keepHills saved in the most recent restart file
Definition: colvarbias_meta.h:195
std::ostringstream hills_traj_os_buf
Cache of the hills trajectory.
Definition: colvarbias_meta.h:297
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:206
int setup_output() override
(Re)initialize the output files (does not write them yet)
Definition: colvarbias_meta.cpp:1742
size_t replica_share_freq() const override
Report the frequency at which this bias needs to communicate with replicas.
Definition: colvarbias_meta.cpp:1019
size_t update_status
Definition: colvarbias_meta.h:286
std::ostream & write_state_data(std::ostream &os) override
Write all mutable data not already written by get_state_params() to a formatted stream.
Definition: colvarbias_meta.cpp:1897
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:268
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:2025
bool b_hills_traj
Write the hill logfile.
Definition: colvarbias_meta.h:110
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:800
std::shared_ptr< colvar_grid_gradient > hills_energy_gradients
Hill forces, cached on a grid.
Definition: colvarbias_meta.h:233
cvm::real hill_width
Definition: colvarbias_meta.h:101
bool replica_state_file_in_sync
Whether a mirror bias has read the latest version of its state file.
Definition: colvarbias_meta.h:282
void rebin_grids_after_restart()
Function called by read_state_data() to execute rebinning (if requested)
Definition: colvarbias_meta.cpp:1481
std::streampos replica_hills_file_pos
Position within replica_hills_file (when reading it)
Definition: colvarbias_meta.h:294
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:392
void recount_hills_off_grid(hill_iter h_first, hill_iter h_last)
Regenerate the hills_off_grid list.
Definition: colvarbias_meta.cpp:983
int init(std::string const &conf) override
Parse config string and (re)initialize.
Definition: colvarbias_meta.cpp:49
std::istream & read_state_data(std::istream &is) override
Read all mutable data not already set by set_state_params() from a formatted stream.
Definition: colvarbias_meta.cpp:1469
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:121
std::string replica_file_name
File containing the paths to the output files from this replica.
Definition: colvarbias_meta.h:246
std::string const hills_traj_file_name() const
Name of the hill logfile.
Definition: colvarbias_meta.cpp:1839
virtual int update_replicas_registry()
Read the existing replicas on registry.
Definition: colvarbias_meta.cpp:1025
std::string replicas_registry
List of replicas (and their output list files)
Definition: colvarbias_meta.h:274
bool dump_fes
Dump the free energy surface (.pmf file) every restartFrequency.
Definition: colvarbias_meta.h:198
hill_iter new_hills_off_grid_begin
Same as new_hills_begin, but for the off-grid ones.
Definition: colvarbias_meta.h:128
std::vector< colvarbias_meta * > replicas
Additional, "mirror" metadynamics biases, to collect info from the other replicas.
Definition: colvarbias_meta.h:265
std::vector< cvm::real > colvar_sigmas
The sigma parameters of the Gaussian hills.
Definition: colvarbias_meta.h:104
Base class for unconstrained thermodynamic-integration FE estimator.
Definition: colvarbias.h:336
Collective variable bias, base class.
Definition: colvarbias.h:23
double real
Defining an abstract real number allows to switch precision.
Definition: colvarmodule.h:150
long long step_number
Use a 64-bit integer to store the step number.
Definition: colvarmodule.h:147
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