Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2012-2019 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_vesselbase_ActionWithVessel_h
23 : #define __PLUMED_vesselbase_ActionWithVessel_h
24 :
25 : #include "core/ActionWithValue.h"
26 : #include "core/ActionAtomistic.h"
27 : #include "tools/Exception.h"
28 : #include "tools/DynamicList.h"
29 : #include "tools/MultiValue.h"
30 : #include <vector>
31 :
32 : namespace PLMD {
33 : class Value;
34 : class Stopwatch;
35 :
36 : namespace vesselbase {
37 :
38 : class Vessel;
39 : class BridgeVessel;
40 : class StoreDataVessel;
41 :
42 : /**
43 : \ingroup MULTIINHERIT
44 : This is used to create PLMD::Action objects that are computed by calculating the same function multiple
45 : times. This is used in PLMD::MultiColvar.
46 : */
47 :
48 : class ActionWithVessel : public virtual Action {
49 : friend class Vessel;
50 : friend class ShortcutVessel;
51 : friend class FunctionVessel;
52 : friend class StoreDataVessel;
53 : friend class BridgeVessel;
54 : friend class ActionWithInputVessel;
55 : friend class OrderingVessel;
56 : private:
57 : /// Do all calculations in serial
58 : bool serial;
59 : /// Lower memory requirements
60 : bool lowmem;
61 : /// Are we skipping the calculation of the derivatives
62 : bool noderiv;
63 : /// This tells plumed that this is used in a bridge
64 : bool actionIsBridged;
65 : /// The maximum number of derivatives we can use before we need to invoke lowmem
66 : unsigned maxderivatives;
67 : /// The tolerance on the accumulators
68 : double tolerance;
69 : /// Tolerance for quantities being put in neighbor lists
70 : double nl_tolerance;
71 : /// Pointers to the functions we are using on each value
72 : std::vector<Vessel*> functions;
73 : /// Tempory storage for forces
74 : std::vector<double> tmpforces;
75 : /// Ths full list of tasks we have to perform
76 : std::vector<unsigned> fullTaskList;
77 : /// The current number of active tasks
78 : unsigned nactive_tasks;
79 : /// The indices of the tasks in the full list of tasks
80 : std::vector<unsigned> indexOfTaskInFullList;
81 : /// The list of currently active tasks
82 : std::vector<unsigned> partialTaskList;
83 : /// The list of atoms involved in derivatives (we keep a copy here to avoid resizing)
84 : std::vector<unsigned> der_list;
85 : /// The buffer that we use (we keep a copy here to avoid resizing)
86 : std::vector<double> buffer;
87 : /// Do we want to output information on the timings of different parts of the calculation
88 : bool timers;
89 : /// The stopwatch that times the different parts of the calculation
90 : Stopwatch& stopwatch;
91 : /// These are used to minmise computational expense in complex functions
92 : bool dertime_can_be_off;
93 : protected:
94 : /// This is also used to minimise computational expense in complex functions
95 : bool dertime;
96 : /// The terms in the series are locked
97 : bool contributorsAreUnlocked;
98 : /// Does the weight have derivatives
99 : bool weightHasDerivatives;
100 : /// This is used for numerical derivatives of bridge variables
101 : unsigned bridgeVariable;
102 : /// A pointer to the object that stores data
103 : StoreDataVessel* mydata;
104 : /// This list is used to update the neighbor list
105 : std::vector<unsigned> taskFlags;
106 : /// Add a vessel to the list of vessels
107 : void addVessel( const std::string& name, const std::string& input, const int numlab=0 );
108 : void addVessel( Vessel* vv );
109 : /// Add a bridging vessel to the list of vessels
110 : BridgeVessel* addBridgingVessel( ActionWithVessel* tome );
111 : /// Complete the setup of this object (this routine must be called after construction of ActionWithValue)
112 : void readVesselKeywords();
113 : /// Turn on the derivatives in the vessel
114 : void needsDerivatives();
115 : /// Return the value of the tolerance
116 : double getTolerance() const ;
117 : /// Return the value for the neighbor list tolerance
118 : double getNLTolerance() const ;
119 : /// Calculate the values of all the vessels
120 : void runAllTasks();
121 : /// Resize all the functions when the number of derivatives change
122 : void resizeFunctions();
123 : /// This loops over all the vessels calculating them and also
124 : /// sets all the element derivatives equal to zero
125 : void calculateAllVessels( const unsigned& taskCode, MultiValue& myvals, MultiValue& bvals, std::vector<double>& buffer, std::vector<unsigned>& der_list );
126 : /// Retrieve the forces from all the vessels (used in apply)
127 : bool getForcesFromVessels( std::vector<double>& forcesToApply );
128 : /// Is the calculation being done in serial
129 : bool serialCalculation() const;
130 : /// Are we using low memory
131 : bool usingLowMem() const ;
132 : /// Set that we are using low memory
133 : void setLowMemOption(const bool& );
134 : /// Deactivate all the tasks in the task list
135 : void deactivateAllTasks();
136 : /// Get the size of the buffer
137 : unsigned getSizeOfBuffer( unsigned& bufsize );
138 : /// Add a task to the full list
139 : void addTaskToList( const unsigned& taskCode );
140 : public:
141 : static void registerKeywords(Keywords& keys);
142 : explicit ActionWithVessel(const ActionOptions&ao);
143 : ~ActionWithVessel();
144 : void lockContributors();
145 : /// Get the number of tasks that are currently active
146 : unsigned getCurrentNumberOfActiveTasks() const ;
147 : /// Check whether or not a particular task is currently active
148 : bool taskIsCurrentlyActive( const unsigned& index ) const ;
149 : /// Are derivatives required for this quantity
150 : bool derivativesAreRequired() const ;
151 : /// Is this action thread safe
152 5206 : virtual bool threadSafe() const { return true; }
153 : /// Finish running all the calculations
154 : virtual void finishComputations( const std::vector<double>& buffer );
155 : /// Are the base quantities periodic
156 : virtual bool isPeriodic()=0;
157 : /// What are the domains of the base quantities
158 : virtual void retrieveDomain( std::string& min, std::string& max);
159 : /// Get the number of derivatives for final calculated quantity
160 : virtual unsigned getNumberOfDerivatives()=0;
161 : /// Get the number of quantities that are calculated during each task
162 : virtual unsigned getNumberOfQuantities() const ;
163 : /// Get the number of vessels
164 : unsigned getNumberOfVessels() const;
165 : /// Get a pointer to the ith vessel
166 : Vessel* getPntrToVessel( const unsigned& i );
167 : /// Do any jobs that are required before the task list is undertaken
168 : virtual void doJobsRequiredBeforeTaskList();
169 : /// Get the full size of the taskList dynamic list
170 : unsigned getFullNumberOfTasks() const ;
171 : /// Get the position of the ith active task in the full list
172 : unsigned getPositionInFullTaskList( const unsigned& ii ) const ;
173 : /// Get the code for the ii th task in the list
174 : unsigned getTaskCode( const unsigned& ii ) const ;
175 : /// Get the ith of the currently active tasks
176 : unsigned getActiveTask( const unsigned& ii ) const ;
177 : /// Calculate one of the functions in the distribution
178 : virtual void performTask( const unsigned&, const unsigned&, MultiValue& ) const=0;
179 : /// Do the task if we have a bridge
180 : virtual void transformBridgedDerivatives( const unsigned& current, MultiValue& invals, MultiValue& outvals ) const;
181 : /// Ensure that data required in other vessels is stored
182 : StoreDataVessel* buildDataStashes( ActionWithVessel* actionThatUses );
183 : /// Apply forces from bridge vessel - this is rarely used - currently only in ActionVolume
184 0 : virtual void applyBridgeForces( const std::vector<double>& bb ) { plumed_error(); }
185 : /// These are overwritten in MultiColvarFunction
186 : // virtual void activateIndexes( const unsigned&, const unsigned&, const std::vector<unsigned>& ){}
187 : /// Return a particular named vessel
188 : Vessel* getVesselWithName( const std::string& mynam );
189 : /// Does the weight have derivatives
190 : bool weightWithDerivatives() const ;
191 : /// Return the position in the current task list
192 : unsigned getPositionInCurrentTaskList( const unsigned& myind ) const ;
193 : /// These normalizes vectors and is used in StoreDataVessel
194 0 : virtual void normalizeVector( std::vector<double>& vals ) const { plumed_error(); }
195 0 : virtual void normalizeVectorDerivatives( MultiValue& myvals ) const { plumed_error(); }
196 : };
197 :
198 : inline
199 : double ActionWithVessel::getTolerance() const {
200 712678 : return tolerance;
201 : }
202 :
203 : inline
204 : double ActionWithVessel::getNLTolerance() const {
205 : return nl_tolerance;
206 : }
207 :
208 : inline
209 : unsigned ActionWithVessel::getNumberOfVessels() const {
210 12530151 : return functions.size();
211 : }
212 :
213 : inline
214 194034 : unsigned ActionWithVessel::getNumberOfQuantities() const {
215 194034 : return 2;
216 : }
217 :
218 : inline
219 : Vessel* ActionWithVessel::getPntrToVessel( const unsigned& i ) {
220 : plumed_dbg_assert( i<functions.size() );
221 3718 : return functions[i];
222 : }
223 :
224 : inline
225 : unsigned ActionWithVessel::getFullNumberOfTasks() const {
226 20230614 : return fullTaskList.size();
227 : }
228 :
229 : inline
230 : unsigned ActionWithVessel::getTaskCode( const unsigned& ii ) const {
231 : plumed_dbg_assert( ii<fullTaskList.size() );
232 7018602 : return fullTaskList[ii];
233 : }
234 :
235 : inline
236 : unsigned ActionWithVessel::getCurrentNumberOfActiveTasks() const {
237 28504 : return nactive_tasks;
238 : }
239 :
240 : inline
241 : unsigned ActionWithVessel::getActiveTask( const unsigned& ii ) const {
242 : plumed_dbg_assert( ii<nactive_tasks );
243 108932 : return partialTaskList[ii];
244 : }
245 :
246 : inline
247 : unsigned ActionWithVessel::getPositionInFullTaskList( const unsigned& ii ) const {
248 : plumed_dbg_assert( ii<nactive_tasks );
249 126939 : return indexOfTaskInFullList[ii];
250 : }
251 :
252 : inline
253 : bool ActionWithVessel::serialCalculation() const {
254 1450 : return serial;
255 : }
256 :
257 : inline
258 : bool ActionWithVessel::usingLowMem() const {
259 : return lowmem;
260 : }
261 :
262 : inline
263 : void ActionWithVessel::setLowMemOption(const bool& l) {
264 33 : lowmem=l;
265 : }
266 :
267 : inline
268 : bool ActionWithVessel::derivativesAreRequired() const {
269 644232 : return !noderiv;
270 : }
271 :
272 : inline
273 : bool ActionWithVessel::weightWithDerivatives() const {
274 13612 : return weightHasDerivatives;
275 : }
276 :
277 : inline
278 218264 : unsigned ActionWithVessel::getPositionInCurrentTaskList( const unsigned& myind ) const {
279 436528 : if( nactive_tasks==fullTaskList.size() ) return myind;
280 :
281 507477438 : for(unsigned i=0; i<nactive_tasks; ++i) {
282 507519400 : if( myind==indexOfTaskInFullList[i] ) return i;
283 : }
284 0 : plumed_merror("requested task is not active");
285 : }
286 :
287 : }
288 : }
289 : #endif
|