Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2015-2020 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.org for more information.
6 :
7 : This file is part of plumed, version 2.
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_adjmat_AdjacencyMatrixBase_h
23 : #define __PLUMED_adjmat_AdjacencyMatrixBase_h
24 :
25 : #include <vector>
26 : #include "core/ActionWithMatrix.h"
27 : #include "tools/LinkCells.h"
28 :
29 : namespace PLMD {
30 : namespace adjmat {
31 :
32 : class AdjacencyMatrixBase : public ActionWithMatrix {
33 : private:
34 : bool nopbc, components, read_one_group;
35 : bool neighbour_list_updated;
36 : LinkCells linkcells, threecells;
37 : std::vector<unsigned> ablocks, threeblocks;
38 : double nl_cut, nl_cut2;
39 : unsigned maxcol;
40 : unsigned nl_stride;
41 : unsigned natoms_per_list;
42 : std::vector<unsigned> nlist;
43 : void setupThirdAtomBlock( const std::vector<AtomNumber>& tc, std::vector<AtomNumber>& t );
44 : protected:
45 : Vector getPosition( const unsigned& indno, MultiValue& myvals ) const ;
46 : void addAtomDerivatives( const unsigned& indno, const Vector& der, MultiValue& myvals ) const ;
47 : void addThirdAtomDerivatives( const unsigned& indno, const Vector& der, MultiValue& myvals ) const ;
48 : void setLinkCellCutoff( const bool& symmetric, const double& lcut, double tcut=-1.0 );
49 : void addBoxDerivatives( const Tensor& vir, MultiValue& myvals ) const ;
50 : public:
51 : static void registerKeywords( Keywords& keys );
52 : explicit AdjacencyMatrixBase(const ActionOptions&);
53 21246780 : bool isAdjacencyMatrix() const override { return true; }
54 : unsigned getNumberOfDerivatives() override ;
55 : unsigned getNumberOfColumns() const override;
56 : void prepare() override;
57 : void getAdditionalTasksRequired( ActionWithVector* action, std::vector<unsigned>& atasks ) override ;
58 : void setupForTask( const unsigned& current, std::vector<unsigned> & indices, MultiValue& myvals ) const override;
59 : // void setupCurrentTaskList() override;
60 : void updateNeighbourList() override ;
61 : unsigned retrieveNeighbours( const unsigned& current, std::vector<unsigned> & indices ) const ;
62 : void performTask( const std::string& controller, const unsigned& index1, const unsigned& index2, MultiValue& myvals ) const override ;
63 : virtual double calculateWeight( const Vector& pos1, const Vector& pos2, const unsigned& natoms, MultiValue& myvals ) const = 0;
64 : void runEndOfRowJobs( const unsigned& ival, const std::vector<unsigned> & indices, MultiValue& myvals ) const override ;
65 : };
66 :
67 : inline
68 : Vector AdjacencyMatrixBase::getPosition( const unsigned& indno, MultiValue& myvals ) const {
69 238830408 : unsigned index = myvals.getIndices()[ indno + myvals.getSplitIndex() ];
70 238830408 : return myvals.getAtomVector()[index];
71 : }
72 :
73 : inline
74 21364880 : void AdjacencyMatrixBase::addAtomDerivatives( const unsigned& indno, const Vector& der, MultiValue& myvals ) const {
75 21364880 : if( doNotCalculateDerivatives() ) return;
76 : plumed_dbg_assert( indno<2 ); unsigned index = myvals.getTaskIndex();
77 12728180 : if( indno==1 ) index = myvals.getSecondTaskIndex();
78 12728180 : unsigned w_index = getConstPntrToComponent(0)->getPositionInStream();
79 12728180 : myvals.addDerivative( w_index, 3*index+0, der[0] );
80 12728180 : myvals.addDerivative( w_index, 3*index+1, der[1] );
81 12728180 : myvals.addDerivative( w_index, 3*index+2, der[2] );
82 : }
83 :
84 : inline
85 131645 : void AdjacencyMatrixBase::addThirdAtomDerivatives( const unsigned& indno, const Vector& der, MultiValue& myvals ) const {
86 131645 : if( doNotCalculateDerivatives() ) return;
87 113597 : unsigned index = myvals.getIndices()[ indno + myvals.getSplitIndex() ];
88 113597 : unsigned w_index = getConstPntrToComponent(0)->getPositionInStream();
89 113597 : myvals.addDerivative( w_index, 3*index+0, der[0] );
90 113597 : myvals.addDerivative( w_index, 3*index+1, der[1] );
91 113597 : myvals.addDerivative( w_index, 3*index+2, der[2] );
92 : }
93 :
94 : inline
95 10682440 : void AdjacencyMatrixBase::addBoxDerivatives( const Tensor& vir, MultiValue& myvals ) const {
96 10682440 : if( doNotCalculateDerivatives() ) return;
97 6364090 : unsigned nbase = 3*getNumberOfAtoms();
98 6364090 : unsigned w_index = getConstPntrToComponent(0)->getPositionInStream();
99 6364090 : myvals.addDerivative( w_index, nbase+0, vir(0,0) );
100 6364090 : myvals.addDerivative( w_index, nbase+1, vir(0,1) );
101 6364090 : myvals.addDerivative( w_index, nbase+2, vir(0,2) );
102 6364090 : myvals.addDerivative( w_index, nbase+3, vir(1,0) );
103 6364090 : myvals.addDerivative( w_index, nbase+4, vir(1,1) );
104 6364090 : myvals.addDerivative( w_index, nbase+5, vir(1,2) );
105 6364090 : myvals.addDerivative( w_index, nbase+6, vir(2,0) );
106 6364090 : myvals.addDerivative( w_index, nbase+7, vir(2,1) );
107 6364090 : myvals.addDerivative( w_index, nbase+8, vir(2,2) );
108 : }
109 :
110 : inline
111 284222 : unsigned AdjacencyMatrixBase::getNumberOfColumns() const {
112 284222 : return maxcol;
113 : }
114 :
115 :
116 : }
117 : }
118 :
119 : #endif
|