All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
SecondaryStructureRMSD.h
Go to the documentation of this file.
1 /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2  Copyright (c) 2013 The plumed team
3  (see the PEOPLE file at the root of the distribution for a list of names)
4 
5  See http://www.plumed-code.org for more information.
6 
7  This file is part of plumed, version 2.0.
8 
9  plumed is free software: you can redistribute it and/or modify
10  it under the terms of the GNU Lesser General Public License as published by
11  the Free Software Foundation, either version 3 of the License, or
12  (at your option) any later version.
13 
14  plumed is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  GNU Lesser General Public License for more details.
18 
19  You should have received a copy of the GNU Lesser General Public License
20  along with plumed. If not, see <http://www.gnu.org/licenses/>.
21 +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
22 #ifndef __PLUMED_secondarystructure_SecondaryStructureRMSD_h
23 #define __PLUMED_secondarystructure_SecondaryStructureRMSD_h
24 
25 #include "core/ActionAtomistic.h"
26 #include "core/ActionWithValue.h"
27 #include "vesselbase/ActionWithVessel.h"
28 #include <vector>
29 
30 namespace PLMD {
31 
32 class SingleDomainRMSD;
33 class DRMSD;
34 class RMSD;
35 
36 namespace secondarystructure {
37 
38 /// Base action for calculating things like AlphRMSD, AntibetaRMSD, etc
39 
41  public ActionAtomistic,
42  public ActionWithValue,
44 {
45 private:
46 /// Are we using pbc
47  bool pbcon;
48 /// Tempory integer to say which refernce configuration is the closest
49  unsigned closest;
50 /// The type of rmsd we are calculating
51  std::string alignType;
52 /// List of all the atoms we require
54 /// The atoms involved in each of the secondary structure segments
55  std::vector< std::vector<unsigned> > colvar_atoms;
56 /// The reference configurations
57  std::vector<RMSD*> secondary_rmsd;
58  std::vector<DRMSD*> secondary_drmsd;
59 /// Stuff for derivatives
60  std::vector< std::vector<Vector> > der;
61  std::vector<Tensor> vir;
62 /// Everything for controlling the updating of neighbor lists
64  unsigned lastUpdate;
66 /// Variables for strands cutoff
68  double s_cutoff;
71 /// Tempory variables for getting positions of atoms and applying forces
72  std::vector<Vector> pos;
73  std::vector<double> forcesToApply;
74 /// Get the index of an atom
75  unsigned getAtomIndex( const unsigned& iatom );
76 protected:
77 /// Get the atoms in the backbone
78  void readBackboneAtoms( const std::vector<std::string>& backnames, std::vector<unsigned>& chain_lengths );
79 /// Add a set of atoms to calculat ethe rmsd from
80  void addColvar( const std::vector<unsigned>& newatoms );
81 /// Set a reference configuration
82  void setSecondaryStructure( std::vector<Vector>& structure, double bondlength, double units );
83 /// Setup a pair of atoms to use for strands cutoff
84  void setAtomsFromStrands( const unsigned& atom1, const unsigned& atom2 );
85 public:
86  static void registerKeywords( Keywords& keys );
88  virtual ~SecondaryStructureRMSD();
90  unsigned getNumberOfDerivatives();
91  void prepare();
92  void calculate();
93  void performTask( const unsigned& j );
94  void clearDerivativesAfterTask( const unsigned& ){}
95  void apply();
96  void mergeDerivatives( const unsigned& , const double& );
97  bool isPeriodic(){ return false; }
98 };
99 
100 inline
102  return colvar_atoms.size();
103 }
104 
105 inline
107  return 3*getNumberOfAtoms()+9;
108 }
109 
110 inline
111 unsigned SecondaryStructureRMSD::getAtomIndex( const unsigned& iatom ){
112  return all_atoms.linkIndex( colvar_atoms[current][iatom] );
113 }
114 
115 }
116 }
117 
118 #endif
bool isPeriodic()
Are the base quantities periodic.
unsigned current
The numerical index of the task we are curently performing.
void setSecondaryStructure(std::vector< Vector > &structure, double bondlength, double units)
Set a reference configuration.
A class for storing a list that changes which members are active as a function of time...
Definition: DynamicList.h:135
std::vector< Vector > pos
Tempory variables for getting positions of atoms and applying forces.
int updateFreq
Everything for controlling the updating of neighbor lists.
void readBackboneAtoms(const std::vector< std::string > &backnames, std::vector< unsigned > &chain_lengths)
Get the atoms in the backbone.
std::string alignType
The type of rmsd we are calculating.
std::vector< std::vector< Vector > > der
Stuff for derivatives.
void addColvar(const std::vector< unsigned > &newatoms)
Add a set of atoms to calculat ethe rmsd from.
unsigned getAtomIndex(const unsigned &iatom)
Get the index of an atom.
Used to create a PLMD::Action that has some scalar or vectorial output that may or may not have some ...
DynamicList< AtomNumber > all_atoms
List of all the atoms we require.
This class holds the keywords and their documentation.
Definition: Keywords.h:36
This class is used to bring the relevant information to the Action constructor.
Definition: Action.h:41
Action used to create objects that access the positions of the atoms from the MD code.
bool align_strands
Variables for strands cutoff.
void performTask(const unsigned &j)
Calculate one of the functions in the distribution.
void mergeDerivatives(const unsigned &, const double &)
unsigned getNumberOfAtoms() const
Get number of available atoms.
void setAtomsFromStrands(const unsigned &atom1, const unsigned &atom2)
Setup a pair of atoms to use for strands cutoff.
std::vector< RMSD * > secondary_rmsd
The reference configurations.
unsigned closest
Tempory integer to say which refernce configuration is the closest.
std::vector< std::vector< unsigned > > colvar_atoms
The atoms involved in each of the secondary structure segments.
void prepare()
Prepare an Action for calculation This can be used by Action if they need some special preparation be...
This is used to create PLMD::Action objects that are computed by calculating the same function multip...
Base action for calculating things like AlphRMSD, AntibetaRMSD, etc.