SSAGES  0.1
A MetaDynamics Package
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Groups Pages
AntiBetaRMSDCV.h
1 
21 #pragma once
22 
23 #include "CollectiveVariable.h"
24 #include "Validator/ObjectRequirement.h"
25 #include "Drivers/DriverException.h"
26 #include "Snapshot.h"
27 #include "schema.h"
28 #include "Utility/ReadBackbone.h"
29 #include <stdexcept>
30 
31 namespace SSAGES
32 {
34 
44  {
45  private:
46 
48  std::vector<int> resids_;
49 
51  //std::vector<int> atomids_;
52  std::vector< std::vector<std::string> > atomids_;
53 
55  std::string refpdb_;
56 
58  std::vector<Vector3> refalpha_;
59 
61  double unitconv_;
62 
64  int mode_;
65 
66  public:
68 
78  AntiBetaRMSDCV(std::vector<int> resids, std::string refpdb, double unitconv, int mode) :
79  resids_(resids), refpdb_(refpdb), unitconv_(unitconv), mode_(mode)
80  {
81  if(resids_.size() != 2 ){
82  throw std::invalid_argument("AntiBetaRMSDCV: Input must designate range of residues with 2 residue numbers.");
83  }
84 
85  resids_.clear();
86 
87  if(resids[0] >= resids[1]){
88  throw std::invalid_argument("AntiBetaRMSDCV: Input must list lower residue index first: please reverse residue range.");
89  } else if(resids[1] - resids[0] < 5){
90  throw std::invalid_argument("AntiBetaRMSDCV: Residue range must span at least 6 residues for secondary structure calculation.");
91  }
92 
93  std::cout << "AntiBetaRMSDCV: Calculating antiparallel beta sheet character from residue " << resids[0] << " to " << resids[1] << "." << std::endl;
94 
95  for(unsigned int i = resids[0]; i <= resids[1]; i++){
96  resids_.push_back(i);
97  }
98  }
99 
101 
104  void Initialize(const Snapshot& snapshot) override
105  {
107 
108  // reference 'ideal' antiparallel beta sheet structure in nanometers
109  // Reference structure taken from values used in Plumed 2.2.0
110  refalpha_.push_back(unitconv_ * Vector3{ .2263, -.3795, .1722}); // N i
111  refalpha_.push_back(unitconv_ * Vector3{ .2493, -.2426, .2263}); // CA
112  refalpha_.push_back(unitconv_ * Vector3{ .3847, -.1838, .1761}); // CB
113  refalpha_.push_back(unitconv_ * Vector3{ .1301, -.1517, .1921}); // C
114  refalpha_.push_back(unitconv_ * Vector3{ .0852, -.1504, .0739}); // O
115  refalpha_.push_back(unitconv_ * Vector3{ .0818, -.0738, .2917}); // N i+1
116  refalpha_.push_back(unitconv_ * Vector3{-.0299, .0243, .2748}); // CA
117  refalpha_.push_back(unitconv_ * Vector3{-.1421, -.0076, .3757}); // CB
118  refalpha_.push_back(unitconv_ * Vector3{ .0273, .1680, .2854}); // C
119  refalpha_.push_back(unitconv_ * Vector3{ .0902, .1993, .3888}); // O
120  refalpha_.push_back(unitconv_ * Vector3{ .0119, .2532, .1813}); // N i+2
121  refalpha_.push_back(unitconv_ * Vector3{ .0683, .3916, .1680}); // CA
122  refalpha_.push_back(unitconv_ * Vector3{ .1580, .3940, .0395}); // CB
123  refalpha_.push_back(unitconv_ * Vector3{-.0394, .5011, .1630}); // C
124  refalpha_.push_back(unitconv_ * Vector3{-.1459, .4814, .0982}); // O
125  refalpha_.push_back(unitconv_ * Vector3{-.2962, .3559, -.1359}); // N j - 2
126  refalpha_.push_back(unitconv_ * Vector3{-.2439, .2526, -.2287}); // CA
127  refalpha_.push_back(unitconv_ * Vector3{-.1189, .3006, -.3087}); // CB
128  refalpha_.push_back(unitconv_ * Vector3{-.2081, .1231, -.1520}); // C
129  refalpha_.push_back(unitconv_ * Vector3{-.1524, .1324, -.0409}); // O
130  refalpha_.push_back(unitconv_ * Vector3{-.2326, .0037, -.2095}); // N j - 1
131  refalpha_.push_back(unitconv_ * Vector3{-.1858, -.1269, -.1554}); // CA
132  refalpha_.push_back(unitconv_ * Vector3{-.3053, -.2199, -.1291}); // CB
133  refalpha_.push_back(unitconv_ * Vector3{-.0869, -.1949, -.2512}); // C
134  refalpha_.push_back(unitconv_ * Vector3{-.1255, -.2070, -.3710}); // O
135  refalpha_.push_back(unitconv_ * Vector3{ .0326, -.2363, -.2072}); // N j
136  refalpha_.push_back(unitconv_ * Vector3{ .1405, -.2992, -.2872}); // CA
137  refalpha_.push_back(unitconv_ * Vector3{ .2699, -.2129, -.2917}); // CB
138  refalpha_.push_back(unitconv_ * Vector3{ .1745, -.4399, -.2330}); // C
139  refalpha_.push_back(unitconv_ * Vector3{ .1899, -.4545, -.1102}); // O
140  }
141 
143 
146  void Evaluate(const Snapshot& snapshot) override
147  {
148  // need atom positions for all atoms in atomids_
149  const auto& pos = snapshot.GetPositions();
150  std::vector<int> groupidx;
151  std::vector<int> double_atomids(atomids_[0].size());
152  std::transform(atomids_[0].begin(), atomids_[0].end(), double_atomids.begin(), [](std::string val) {return std::stod(val);});
153  snapshot.GetLocalIndices(double_atomids, &groupidx);
154 
155  unsigned int resgroups = resids_.size() - 2;
156  double rmsd, dist_norm, dxgrouprmsd;
157  Vector3 dist_xyz, dist_ref;
158  std::vector<Vector3> refxyz;
159  std::vector< std::vector< Vector3 > > deriv(30, std::vector<Vector3>(30, Vector3{0,0,0}));
160 
161  std::fill(grad_.begin(), grad_.end(), Vector3{0,0,0});
162  grad_.resize(snapshot.GetNumAtoms(), Vector3{0,0,0});
163  val_ = 0.0;
164 
165  for(size_t i = 0; i < resgroups - 3; i++){
166  for(size_t j = i + 3; j < resgroups; j++){
167  // mode: 0 for all, 1 for inter, 2 for intra
168  if((mode_ == 0) || (mode_ == 1 && atomids_[1][5 * j] != atomids_[1][5 * i]) || (mode_ == 2 && atomids_[1][5 * j] == atomids_[1][5 * i])){
169  rmsd = 0.0;
170  std::fill(refxyz.begin(), refxyz.end(), Vector3{0,0,0});
171  refxyz.resize(30, Vector3{0,0,0});
172  for(size_t k = 0; k < 15; k++){
173  refxyz[k] = pos[groupidx[5 * i + k]];
174  }
175  for(size_t k = 0; k < 15; k++){
176  refxyz[k + 15] = pos[groupidx[5 * j + k]];
177  }
178 
179  // sum over all pairs to calculate CV
180  for(size_t k = 0; k < 29; k++){
181  for(size_t h = k + 1; h < 30; h++){
182  dist_xyz = refxyz[k] - refxyz[h];
183  dist_ref = refalpha_[k] - refalpha_[h];
184  dist_norm = dist_xyz.norm() - dist_ref.norm();
185  rmsd += dist_norm * dist_norm;
186  deriv[k][h] = dist_xyz * dist_norm / dist_xyz.norm();
187  }
188  }
189 
190  rmsd = pow(rmsd / 435, 0.5) / 0.1;
191  val_ += (1 - pow(rmsd, 8.0)) / (1 - pow(rmsd, 12.0));
192 
193  dxgrouprmsd = pow(rmsd, 11.0) + pow(rmsd, 7.0);
194  dxgrouprmsd /= pow(rmsd, 8.0) + pow(rmsd, 4.0) + 1;
195  dxgrouprmsd /= pow(rmsd, 8.0) + pow(rmsd, 4.0) + 1;
196  dxgrouprmsd *= -40. / 435;
197 
198  for(size_t k = 0; k < 15; k++){
199  for(size_t h = k + 1; h < 15; h++){
200  grad_[groupidx[5 * i + k]] += dxgrouprmsd * deriv[k][h];
201  grad_[groupidx[5 * i + h]] -= dxgrouprmsd * deriv[k][h];
202  }
203  for(size_t h = 0; h < 15; h++){
204  grad_[groupidx[5 * i + k]] += dxgrouprmsd * deriv[k][h+15];
205  grad_[groupidx[5 * j + h]] -= dxgrouprmsd * deriv[k][h+15];
206  }
207  }
208  for(size_t k = 0; k < 14; k++){
209  for(size_t h = k + 1; h < 15; h++){
210  grad_[groupidx[5 * j + k]] += dxgrouprmsd * deriv[k+15][h+15];
211  grad_[groupidx[5 * j + h]] -= dxgrouprmsd * deriv[k+15][h+15];
212  }
213  }
214  }
215  }
216  }
217  }
218 
220  static AntiBetaRMSDCV* Build(const Json::Value& json, const std::string& path)
221  {
222  Json::ObjectRequirement validator;
223  Json::Value schema;
224  Json::Reader reader;
225 
226  reader.parse(JsonSchema::AntiBetaRMSDCV, schema);
227  validator.Parse(schema, path);
228 
229  //Validate inputs
230  validator.Validate(json, path);
231  if(validator.HasErrors())
232  throw BuildException(validator.GetErrors());
233 
234  std::vector<int> resids;
235  for(auto& s : json["residue_ids"])
236  resids.push_back(s.asInt());
237  auto reference = json["reference"].asString();
238 
239  double unitconv = json.get("length_unit", 1).asDouble();
240 
241  int mode = json.get("mode", 0).asInt();
242 
243  return new AntiBetaRMSDCV(resids, reference, unitconv, mode);
244  }
245  };
246 }
bool HasErrors()
Check if errors have occured.
Definition: Requirement.h:86
double unitconv_
Length unit conversion: convert 1 nm to your internal MD units (ex. if using angstroms use 10) ...
static AntiBetaRMSDCV * Build(const Json::Value &json, const std::string &path)
Set up collective variable.
unsigned GetNumAtoms() const
Get number of atoms in this snapshot.
Definition: Snapshot.h:200
std::vector< Vector3 > refalpha_
Coordinates for reference structure.
Class containing a snapshot of the current simulation in time.
Definition: Snapshot.h:43
virtual void Parse(Value json, const std::string &path) override
Parse JSON value to generate Requirement(s).
std::vector< Vector3 > grad_
Gradient vector dCv/dxi.
void GetLocalIndices(const Label &ids, Label *indices) const
Definition: Snapshot.h:573
void Initialize(const Snapshot &snapshot) override
Initialize necessary variables.
Exception to be thrown when building the Driver fails.
std::vector< std::string > GetErrors()
Get list of error messages.
Definition: Requirement.h:92
Collective variable to measure antiparallel beta secondary structure.
std::vector< int > resids_
Residue IDs for secondary structure calculation.
Requirements on an object.
std::vector< std::vector< std::string > > atomids_
Atom IDs for secondary structure calculation: backbone of resids_.
std::string refpdb_
Name of pdb reference for system.
Eigen::Vector3d Vector3
Three-dimensional vector.
Definition: types.h:33
double val_
Current value of CV.
Abstract class for a collective variable.
void Evaluate(const Snapshot &snapshot) override
Evaluate the CV.
static std::vector< std::vector< std::string > > GetPdbBackbone(std::string refpdb, std::vector< int > resids)
Read protein backbone atoms from reference PDB file.
Definition: ReadBackbone.h:57
AntiBetaRMSDCV(std::vector< int > resids, std::string refpdb, double unitconv, int mode)
Constructor.
const std::vector< Vector3 > & GetPositions() const
Access the particle positions.
Definition: Snapshot.h:323
int mode_
Specify whether to calculate beta sheet character in intra or inter mode: 0 for either, 1 for inter, 2 for intra.
virtual void Validate(const Value &json, const std::string &path) override
Validate JSON value.