File: colvaratoms.h

package info (click to toggle)
lammps 0~20181211.gitad1b1897d%2Bdfsg1-2
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 318,860 kB
  • sloc: cpp: 729,569; python: 40,508; xml: 14,919; fortran: 12,142; ansic: 7,454; sh: 5,565; perl: 4,105; f90: 2,700; makefile: 2,117; objc: 238; lisp: 163; tcl: 61; csh: 16; awk: 14
file content (509 lines) | stat: -rw-r--r-- 15,827 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
// -*- c++ -*-

// This file is part of the Collective Variables module (Colvars).
// The original version of Colvars and its updates are located at:
// https://github.com/colvars/colvars
// Please update all Colvars source files before making any changes.
// If you wish to distribute your changes, please submit them to the
// Colvars repository at GitHub.

#ifndef COLVARATOMS_H
#define COLVARATOMS_H

#include "colvarmodule.h"
#include "colvarproxy.h"
#include "colvarparse.h"
#include "colvardeps.h"


/// \brief Stores numeric id, mass and all mutable data for an atom,
/// mostly used by a \link cvc \endlink
///
/// This class may be used to keep atomic data such as id, mass,
/// position and collective variable derivatives) altogether.
/// There may be multiple instances with identical
/// numeric id, all acting independently: forces communicated through
/// these instances will be summed together.

class colvarmodule::atom {

protected:

  /// Index in the colvarproxy arrays (\b NOT in the global topology!)
  int index;

public:

  /// Identifier for the MD program (0-based)
  int             id;

  /// Mass
  cvm::real       mass;

  /// Charge
  cvm::real       charge;

  /// \brief Current position (copied from the program, can be
  /// modified if necessary)
  cvm::atom_pos   pos;

  /// \brief Current velocity (copied from the program, can be
  /// modified if necessary)
  cvm::rvector    vel;

  /// \brief System force at the previous step (copied from the
  /// program, can be modified if necessary)
  cvm::rvector    total_force;

  /// \brief Gradient of a scalar collective variable with respect
  /// to this atom
  ///
  /// This can only handle a scalar collective variable (i.e. when
  /// the \link colvarvalue::real_value \endlink member is used
  /// from the \link colvarvalue \endlink class), which is also the
  /// most frequent case. For more complex types of \link
  /// colvarvalue \endlink objects, atomic gradients should be
  /// defined within the specific \link cvc \endlink
  /// implementation
  cvm::rvector   grad;

  /// \brief Default constructor (sets index and id both to -1)
  atom();

  /// \brief Initialize an atom for collective variable calculation
  /// and get its internal identifier \param atom_number Atom index in
  /// the system topology (1-based)
  atom(int atom_number);

  /// \brief Initialize an atom for collective variable calculation
  /// and get its internal identifier \param residue Residue number
  /// \param atom_name Name of the atom in the residue \param
  /// segment_id For PSF topologies, the segment identifier; for other
  /// type of topologies, may not be required
  atom(cvm::residue_id const &residue,
       std::string const     &atom_name,
       std::string const     &segment_id);

  /// Copy constructor
  atom(atom const &a);

  /// Destructor
  ~atom();

  /// Set mutable data (everything except id and mass) to zero
  inline void reset_data()
  {
    pos = cvm::atom_pos(0.0);
    vel = grad = total_force = cvm::rvector(0.0);
  }

  /// Get the latest value of the mass
  inline void update_mass()
  {
    mass = (cvm::proxy)->get_atom_mass(index);
  }

  /// Get the latest value of the charge
  inline void update_charge()
  {
    charge = (cvm::proxy)->get_atom_charge(index);
  }

  /// Get the current position
  inline void read_position()
  {
    pos = (cvm::proxy)->get_atom_position(index);
  }

  /// Get the current velocity
  inline void read_velocity()
  {
    vel = (cvm::proxy)->get_atom_velocity(index);
  }

  /// Get the total force
  inline void read_total_force()
  {
    total_force = (cvm::proxy)->get_atom_total_force(index);
  }

  /// \brief Apply a force to the atom
  ///
  /// Note: the force is not applied instantly, but will be used later
  /// by the MD integrator (the colvars module does not integrate
  /// equations of motion.
  ///
  /// Multiple calls to this function by either the same
  /// \link atom \endlink object or different objects with identical
  /// \link id \endlink will all be added together.
  inline void apply_force(cvm::rvector const &new_force) const
  {
    (cvm::proxy)->apply_atom_force(index, new_force);
  }
};



/// \brief Group of \link atom \endlink objects, mostly used by a
/// \link cvc \endlink object to gather all atomic data
class colvarmodule::atom_group
  : public colvarparse, public colvardeps
{
public:


  /// \brief Default constructor
  atom_group();

  /// \brief Create a group object, assign a name to it
  atom_group(char const *key);

  /// \brief Initialize the group after a (temporary) vector of atoms
  atom_group(std::vector<cvm::atom> const &atoms_in);

  /// \brief Destructor
  ~atom_group();

  /// \brief Optional name to reuse properties of this in other groups
  std::string name;

  /// \brief Keyword used to define the group
  // TODO Make this field part of the data structures that link a group to a CVC
  std::string key;

  /// \brief Set default values for common flags
  int init();

  /// \brief Update data required to calculate cvc's
  int setup();

  /// \brief Initialize the group by looking up its configuration
  /// string in conf and parsing it
  int parse(std::string const &conf);

  int add_atom_numbers(std::string const &numbers_conf);
  int add_atoms_of_group(atom_group const * ag);
  int add_index_group(std::string const &index_group_name);
  int add_atom_numbers_range(std::string const &range_conf);
  int add_atom_name_residue_range(std::string const &psf_segid,
                                  std::string const &range_conf);
  int parse_fitting_options(std::string const &group_conf);

  /// \brief Add an atom object to this group
  int add_atom(cvm::atom const &a);

  /// \brief Add an atom ID to this group (the actual atomicdata will be not be handled by the group)
  int add_atom_id(int aid);

  /// \brief Remove an atom object from this group
  int remove_atom(cvm::atom_iter ai);

  /// \brief Re-initialize the total mass of a group.
  /// This is needed in case the hosting MD code has an option to
  /// change atom masses after their initialization.
  void reset_mass(std::string &name, int i, int j);

  /// \brief Implementation of the feature list for atom group
  static std::vector<feature *> ag_features;

  /// \brief Implementation of the feature list accessor for atom group
  virtual const std::vector<feature *> &features()
  {
    return ag_features;
  }
  virtual std::vector<feature *> &modify_features()
  {
    return ag_features;
  }
  static void delete_features() {
    for (size_t i=0; i < ag_features.size(); i++) {
      delete ag_features[i];
    }
    ag_features.clear();
  }

protected:

  /// \brief Array of atom objects
  std::vector<cvm::atom> atoms;

  /// \brief Internal atom IDs for host code
  std::vector<int> atoms_ids;

  /// Sorted list of internal atom IDs (populated on-demand by
  /// create_sorted_ids); used to read coordinate files
  std::vector<int> sorted_atoms_ids;

  /// Map entries of sorted_atoms_ids onto the original positions in the group
  std::vector<int> sorted_atoms_ids_map;

  /// \brief Dummy atom position
  cvm::atom_pos dummy_atom_pos;

  /// \brief Index in the colvarproxy arrays (if the group is scalable)
  int index;

public:

  inline cvm::atom & operator [] (size_t const i)
  {
    return atoms[i];
  }

  inline cvm::atom const & operator [] (size_t const i) const
  {
    return atoms[i];
  }

  inline cvm::atom_iter begin()
  {
    return atoms.begin();
  }

  inline cvm::atom_const_iter begin() const
  {
    return atoms.begin();
  }

  inline cvm::atom_iter end()
  {
    return atoms.end();
  }

  inline cvm::atom_const_iter end() const
  {
    return atoms.end();
  }

  inline size_t size() const
  {
    return atoms.size();
  }

  /// \brief If this option is on, this group merely acts as a wrapper
  /// for a fixed position; any calls to atoms within or to
  /// functions that return disaggregated data will fail
  bool b_dummy;

  /// Internal atom IDs (populated during initialization)
  inline std::vector<int> const &ids() const
  {
    return atoms_ids;
  }

  std::string const print_atom_ids() const;

  /// Allocates and populates sorted_ids and sorted_ids_map
  int create_sorted_ids();

  /// Sorted internal atom IDs (populated on-demand by create_sorted_ids);
  /// used to read coordinate files
  inline std::vector<int> const &sorted_ids() const
  {
    return sorted_atoms_ids;
  }

  /// Map entries of sorted_atoms_ids onto the original positions in the group
  inline std::vector<int> const &sorted_ids_map() const
  {
    return sorted_atoms_ids_map;
  }

  /// Detect whether two groups share atoms
  /// If yes, returns 1-based number of a common atom; else, returns 0
  static int overlap(const atom_group &g1, const atom_group &g2);

  /// \brief When updating atomic coordinates, translate them to align with the
  /// center of mass of the reference coordinates
  bool b_center;

  /// \brief When updating atom coordinates (and after
  /// centering them if b_center is set), rotate the group to
  /// align with the reference coordinates.
  ///
  /// Note: gradients will be calculated in the rotated frame: when
  /// forces will be applied, they will rotated back to the original
  /// frame
  bool b_rotate;
  /// The rotation calculated automatically if b_rotate is defined
  cvm::rotation rot;

  /// \brief Indicates that the user has explicitly set centerReference or
  /// rotateReference, and the corresponding reference:
  /// cvc's (eg rmsd, eigenvector) will not override the user's choice
  bool b_user_defined_fit;

  /// \brief use reference coordinates for b_center or b_rotate
  std::vector<cvm::atom_pos> ref_pos;

  /// \brief Center of geometry of the reference coordinates; regardless
  /// of whether b_center is true, ref_pos is centered to zero at
  /// initialization, and ref_pos_cog serves to center the positions
  cvm::atom_pos              ref_pos_cog;

  /// \brief If b_center or b_rotate is true, use this group to
  /// define the transformation (default: this group itself)
  atom_group                *fitting_group;

  /// Total mass of the atom group
  cvm::real total_mass;
  void update_total_mass();

  /// Total charge of the atom group
  cvm::real total_charge;
  void update_total_charge();

  /// \brief Don't apply any force on this group (use its coordinates
  /// only to calculate a colvar)
  bool        noforce;

  /// \brief Get the current positions
  void read_positions();

  /// \brief (Re)calculate the optimal roto-translation
  void calc_apply_roto_translation();

  /// \brief Save aside the center of geometry of the reference positions,
  /// then subtract it from them
  ///
  /// In this way it will be possible to use ref_pos also for the
  /// rotational fit.
  /// This is called either by atom_group::parse or by CVCs that assign
  /// reference positions (eg. RMSD, eigenvector).
  void center_ref_pos();

  /// \brief Move all positions
  void apply_translation(cvm::rvector const &t);

  /// \brief Get the current velocities; this must be called always
  /// *after* read_positions(); if b_rotate is defined, the same
  /// rotation applied to the coordinates will be used
  void read_velocities();

  /// \brief Get the current total_forces; this must be called always
  /// *after* read_positions(); if b_rotate is defined, the same
  /// rotation applied to the coordinates will be used
  void read_total_forces();

  /// Call reset_data() for each atom
  inline void reset_atoms_data()
  {
    for (cvm::atom_iter ai = atoms.begin(); ai != atoms.end(); ai++)
      ai->reset_data();
    if (fitting_group)
      fitting_group->reset_atoms_data();
  }

  /// \brief Recompute all mutable quantities that are required to compute CVCs
  int calc_required_properties();

  /// \brief Return a copy of the current atom positions
  std::vector<cvm::atom_pos> positions() const;

  /// \brief Calculate the center of geometry of the atomic positions, assuming
  /// that they are already pbc-wrapped
  int calc_center_of_geometry();

private:

  /// \brief Center of geometry
  cvm::atom_pos cog;

  /// \brief Center of geometry before any fitting
  cvm::atom_pos cog_orig;

public:

  /// \brief Return the center of geometry of the atomic positions
  inline cvm::atom_pos center_of_geometry() const
  {
    return cog;
  }

  /// \brief Calculate the center of mass of the atomic positions, assuming that
  /// they are already pbc-wrapped
  int calc_center_of_mass();
private:
  /// \brief Center of mass
  cvm::atom_pos com;
  /// \brief The derivative of a scalar variable with respect to the COM
  // TODO for scalable calculations of more complex variables (e.g. rotation),
  // use a colvarvalue of vectors to hold the entire derivative
  cvm::rvector scalar_com_gradient;
public:
  /// \brief Return the center of mass of the atomic positions
  inline cvm::atom_pos center_of_mass() const
  {
    return com;
  }

  /// \brief Return a copy of the current atom positions, shifted by a constant vector
  std::vector<cvm::atom_pos> positions_shifted(cvm::rvector const &shift) const;

  /// \brief Return a copy of the current atom velocities
  std::vector<cvm::rvector> velocities() const;

  ///\brief Calculate the dipole of the atom group around the specified center
  int calc_dipole(cvm::atom_pos const &com);
private:
  cvm::rvector dip;
public:
  ///\brief Return the (previously calculated) dipole of the atom group
  inline cvm::rvector dipole() const
  {
    return dip;
  }

  /// \brief Return a copy of the total forces
  std::vector<cvm::rvector> total_forces() const;

  /// \brief Return a copy of the aggregated total force on the group
  cvm::rvector total_force() const;


  /// \brief Shorthand: save the specified gradient on each atom,
  /// weighting with the atom mass (mostly used in combination with
  /// \link center_of_mass() \endlink)
  void set_weighted_gradient(cvm::rvector const &grad);

  /// \brief Calculate the derivatives of the fitting transformation
  void calc_fit_gradients();

  /// \brief Derivatives of the fitting transformation
  std::vector<cvm::atom_pos> fit_gradients;

  /// \brief Used by a (scalar) colvar to apply its force on its \link
  /// atom_group \endlink members
  ///
  /// The (scalar) force is multiplied by the colvar gradient for each
  /// atom; this should be used when a colvar with scalar \link
  /// colvarvalue \endlink type is used (this is the most frequent
  /// case: for colvars with a non-scalar type, the most convenient
  /// solution is to sum together the Cartesian forces from all the
  /// colvar components, and use apply_force() or apply_forces()).  If
  /// the group is being rotated to a reference frame (e.g. to express
  /// the colvar independently from the solute rotation), the
  /// gradients are temporarily rotated to the original frame.
  void apply_colvar_force(cvm::real const &force);

  /// \brief Apply a force "to the center of mass", i.e. the force is
  /// distributed on each atom according to its mass
  ///
  /// If the group is being rotated to a reference frame (e.g. to
  /// express the colvar independently from the solute rotation), the
  /// force is rotated back to the original frame.  Colvar gradients
  /// are not used, either because they were not defined (e.g because
  /// the colvar has not a scalar value) or the biases require to
  /// micromanage the force.
  /// This function will be phased out eventually, in favor of
  /// apply_colvar_force() once that is implemented for non-scalar values
  void apply_force(cvm::rvector const &force);

  /// Implements possible actions to be carried out
  /// when a given feature is enabled
  /// This overloads the base function in colvardeps
  void do_feature_side_effects(int id);
};


#endif