All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
MultiColvarFunction.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_multicolvar_MultiColvarFunction_h
23 #define __PLUMED_multicolvar_MultiColvarFunction_h
24 
25 #include "MultiColvarBase.h"
27 
28 namespace PLMD {
29 namespace multicolvar {
30 
32 private:
33 /// The multicolvar from which we construct these quantities
35 /// The central atom positions
37 protected:
38 /// Get the number of functions in the multicolvar we are operating on
39  unsigned getNumberOfBaseFunctions() const;
40 /// Return a pointer to the multicolvar we are using as a base function
42 /// Finish off the setup of the VectorFunction
43  void completeSetup();
44 /// Get the position of one of the central atoms
45  Vector getPositionOfCentralAtom(const unsigned&) const;
46 /// Add derivatives of value wrt to an atomic position
47  void addCentralAtomsDerivatives( const unsigned& , const Vector& );
48 /// Add derivatives of weight wrt to an atomic position
49  void addCentralAtomsDerivativesOfWeight( const unsigned& , const Vector& );
50 /// We use the distance from mycolv as otherwise numerical derivatives dont work
51  Vector getSeparation( const Vector& vec1, const Vector& vec2 ) const ;
52 /// Add derivative wrt to the position of the central atom
53  void addDerivativeOfCentralAtomPos( const unsigned& iatom, const Tensor& der );
54 public:
56  static void registerKeywords( Keywords& keys );
57 /// Used to make sure we are calculating everything during neighbor list update step
58  void unlockContributors();
59  void lockContributors();
60 /// Active element in atoms_with_derivatives
61  void atomHasDerivative( const unsigned& iatom );
62 /// Resize the dynamic arrays
63  void resizeDynamicArrays();
64 /// Regular calculate
65  void calculate();
66 /// Calculate the numerical derivatives for this action
68 /// Do the calculation
69  double doCalculation( const unsigned& j );
70 /// Actually compute the colvar
71  virtual double compute( const unsigned& j )=0;
72 /// Calculate the position of the central atom
74 /// Get the position of the central atom
75  virtual Vector getCentralAtom()=0;
76 };
77 
78 inline
80  return mycolv->colvar_atoms.size();
81 }
82 
83 inline
85  return mycolv;
86 }
87 
88 inline
89 Vector MultiColvarFunction::getPositionOfCentralAtom( const unsigned& iatom ) const {
90  plumed_dbg_assert( iatom<colvar_atoms[current].getNumberActive() );
91  return catoms->getPosition( colvar_atoms[current][iatom] );
92 }
93 
94 inline
95 void MultiColvarFunction::addCentralAtomsDerivatives( const unsigned& iatom, const Vector& der ){
96  plumed_dbg_assert( iatom<colvar_atoms[current].getNumberActive() );
97  catoms->addAtomsDerivatives( colvar_atoms[current][iatom], der, this );
98 }
99 
100 inline
101 void MultiColvarFunction::addCentralAtomsDerivativesOfWeight( const unsigned& iatom, const Vector& der ){
102  plumed_dbg_assert( iatom<colvar_atoms[current].getNumberActive() );
103  catoms->addAtomsDerivativeOfWeight( colvar_atoms[current][iatom], der, this );
104 }
105 
106 inline
107 void MultiColvarFunction::atomHasDerivative( const unsigned& iatom ){
109 }
110 
111 inline
112 void MultiColvarFunction::addDerivativeOfCentralAtomPos( const unsigned& iatom, const Tensor& der ){
114 }
115 
116 
117 }
118 }
119 #endif
void addCentralAtomsDerivativesOfWeight(const unsigned &, const Vector &)
Add derivatives of weight wrt to an atomic position.
Vector getSeparation(const Vector &vec1, const Vector &vec2) const
We use the distance from mycolv as otherwise numerical derivatives dont work.
Vector calculateCentralAtomPosition()
Calculate the position of the central atom.
unsigned current
The numerical index of the task we are curently performing.
void addAtomsDerivativeOfWeight(const unsigned &iatom, const Vector &df, MultiColvarFunction *funcout) const
Add derivatives of the weight wrt to the central atom position.
Class implementing fixed size matrices of doubles.
Definition: Tensor.h:70
Class implementing fixed size vectors of doubles.
Definition: Vector.h:74
void addAtomsDerivatives(const unsigned &iatom, const Vector &df, MultiColvarFunction *funcout) const
Add derivatives to central atom position.
void unlockContributors()
Used to make sure we are calculating everything during neighbor list update step. ...
MultiColvarBase * getPntrToMultiColvar()
Return a pointer to the multicolvar we are using as a base function.
multicolvar::StoreCentralAtomsVessel * catoms
The central atom positions.
void activate(const unsigned ii)
Make something active.
Definition: DynamicList.h:263
void addDerivativeOfCentralAtomPos(const unsigned &iatom, const Tensor &df, MultiColvarFunction *funcout) const
Add derivative to the central atom position.
virtual Vector getCentralAtom()=0
Get the position of the central atom.
Used to create a PLMD::Action that has some scalar or vectorial output that may or may not have some ...
void addDerivativeOfCentralAtomPos(const unsigned &iatom, const Tensor &der)
Add derivative wrt to the position of the central atom.
This class holds the keywords and their documentation.
Definition: Keywords.h:36
virtual double compute(const unsigned &j)=0
Actually compute the colvar.
std::vector< DynamicList< unsigned > > colvar_atoms
The lists of the atoms involved in each of the individual colvars note these refer to the atoms in al...
This class is used to bring the relevant information to the Action constructor.
Definition: Action.h:41
void addCentralAtomsDerivatives(const unsigned &, const Vector &)
Add derivatives of value wrt to an atomic position.
void completeSetup()
Finish off the setup of the VectorFunction.
void resizeDynamicArrays()
Resize the dynamic arrays.
Vector getPosition(const unsigned &) const
Get the orientation of the ith vector.
multicolvar::MultiColvarBase * mycolv
The multicolvar from which we construct these quantities.
void calculateNumericalDerivatives(ActionWithValue *a=NULL)
Calculate the numerical derivatives for this action.
double doCalculation(const unsigned &j)
Do the calculation.
unsigned getNumberOfBaseFunctions() const
Get the number of functions in the multicolvar we are operating on.
DynamicList< unsigned > atoms_with_derivatives
A dynamic list containing those atoms with derivatives.
static void registerKeywords(Keywords &keys)
void const char const char int double * a
Definition: Matrix.h:42
Vector getPositionOfCentralAtom(const unsigned &) const
Get the position of one of the central atoms.
void atomHasDerivative(const unsigned &iatom)
Active element in atoms_with_derivatives.