Merge release-4-6 (commit 'Ic142a690')
[alexxy/gromacs.git] / src / gromacs / analysisdata / datastorage.cpp
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11  * Copyright (c) 2001-2009, The GROMACS development team,
12  * check out http://www.gromacs.org for more information.
13
14  * This program is free software; you can redistribute it and/or
15  * modify it under the terms of the GNU General Public License
16  * as published by the Free Software Foundation; either version 2
17  * of the License, or (at your option) any later version.
18  *
19  * If you want to redistribute modifications, please consider that
20  * scientific software is very special. Version control is crucial -
21  * bugs must be traceable. We will be happy to consider code for
22  * inclusion in the official distribution, but derived work must not
23  * be called official GROMACS. Details are found in the README & COPYING
24  * files - if they are missing, get the official version at www.gromacs.org.
25  *
26  * To help us fund GROMACS development, we humbly ask that you cite
27  * the papers on the package - you can find them in the top README file.
28  *
29  * For more info, check our website at http://www.gromacs.org
30  */
31 /*! \internal \file
32  * \brief
33  * Implements classes in datastorage.h and paralleloptions.h.
34  *
35  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36  * \ingroup module_analysisdata
37  */
38 #include "datastorage.h"
39
40 #include <limits>
41
42 #include "gromacs/analysisdata/abstractdata.h"
43 #include "gromacs/analysisdata/dataframe.h"
44 #include "gromacs/analysisdata/paralleloptions.h"
45 #include "gromacs/utility/exceptions.h"
46 #include "gromacs/utility/gmxassert.h"
47
48 #include "datastorage-impl.h"
49
50 namespace gmx
51 {
52
53 /********************************************************************
54  * AnalysisDataParallelOptions
55  */
56
57 AnalysisDataParallelOptions::AnalysisDataParallelOptions()
58     : parallelizationFactor_(1)
59 {
60 }
61
62
63 AnalysisDataParallelOptions::AnalysisDataParallelOptions(int parallelizationFactor)
64     : parallelizationFactor_(parallelizationFactor)
65 {
66     GMX_RELEASE_ASSERT(parallelizationFactor >= 1,
67                        "Invalid parallelization factor");
68 }
69
70
71 /********************************************************************
72  * AnalysisDataStorageFrame
73  */
74
75 AnalysisDataStorageFrame::AnalysisDataStorageFrame(AnalysisDataStorage *storage,
76                                                    int columnCount, int index)
77     : storage_(*storage), header_(index, 0.0, 0.0), values_(columnCount)
78 {
79 }
80
81
82 AnalysisDataStorageFrame::~AnalysisDataStorageFrame()
83 {
84 }
85
86
87 AnalysisDataPointSetRef
88 AnalysisDataStorageFrame::currentPoints() const
89 {
90     std::vector<AnalysisDataValue>::const_iterator begin = values_.begin();
91     std::vector<AnalysisDataValue>::const_iterator end = values_.end();
92     while (begin != end && !begin->isSet())
93     {
94         ++begin;
95     }
96     while (end != begin && !(end-1)->isSet())
97     {
98         --end;
99     }
100     int firstColumn = (begin != end) ? begin - values_.begin() : 0;
101     return AnalysisDataPointSetRef(header_, firstColumn,
102                 AnalysisDataValuesRef(begin, end));
103 }
104
105
106 void
107 AnalysisDataStorageFrame::clearValues()
108 {
109     std::vector<AnalysisDataValue>::iterator i;
110     for (i = values_.begin(); i != values_.end(); ++i)
111     {
112         i->clear();
113     }
114 }
115
116
117 void
118 AnalysisDataStorageFrame::finishPointSet()
119 {
120     storage_.impl_->notifyPointSet(currentPoints());
121     clearValues();
122 }
123
124
125 /********************************************************************
126  * AnalysisDataStorage::Impl
127  */
128
129 AnalysisDataStorage::Impl::Impl()
130     : data_(NULL), bMultipoint_(false),
131       storageLimit_(0), pendingLimit_(1), firstFrameLocation_(0), nextIndex_(0)
132 {
133 }
134
135
136 AnalysisDataStorage::Impl::~Impl()
137 {
138 }
139
140
141 int
142 AnalysisDataStorage::Impl::columnCount() const
143 {
144     GMX_ASSERT(data_ != NULL, "columnCount() called too early");
145     return data_->columnCount();
146 }
147
148
149 bool
150 AnalysisDataStorage::Impl::isMultipoint() const
151 {
152     return bMultipoint_;
153 }
154
155
156 int
157 AnalysisDataStorage::Impl::firstStoredIndex() const
158 {
159     return frames_[firstFrameLocation_].frame->frameIndex();
160 }
161
162
163 int
164 AnalysisDataStorage::Impl::computeStorageLocation(int index) const
165 {
166     if (index < firstStoredIndex() || index >= nextIndex_)
167     {
168         return -1;
169     }
170     return index % frames_.size();
171 }
172
173
174 size_t
175 AnalysisDataStorage::Impl::endStorageLocation() const
176 {
177     if (storeAll())
178     {
179         return frames_.size();
180     }
181     if (frames_[0].frame->frameIndex() == 0 || firstFrameLocation_ == 0)
182     {
183         return frames_.size() - 1;
184     }
185     return firstFrameLocation_ - 1;
186 }
187
188
189 void
190 AnalysisDataStorage::Impl::extendBuffer(AnalysisDataStorage *storage,
191                                         size_t newSize)
192 {
193     frames_.reserve(newSize);
194     while (frames_.size() < newSize)
195     {
196         frames_.push_back(StoredFrame(
197             new AnalysisDataStorageFrame(storage, columnCount(), nextIndex_)));
198         ++nextIndex_;
199     }
200     // The unused frame should not be included in the count.
201     if (!storeAll())
202     {
203         --nextIndex_;
204     }
205 }
206
207
208 void
209 AnalysisDataStorage::Impl::rotateBuffer()
210 {
211     GMX_ASSERT(!storeAll(),
212                "No need to rotate internal buffer if everything is stored");
213     size_t prevFirst = firstFrameLocation_;
214     size_t nextFirst = prevFirst + 1;
215     if (nextFirst == frames_.size())
216     {
217         nextFirst = 0;
218     }
219     firstFrameLocation_ = nextFirst;
220     StoredFrame &prevFrame = frames_[prevFirst];
221     prevFrame.status = StoredFrame::eMissing;
222     prevFrame.frame->header_ = AnalysisDataFrameHeader(nextIndex_ + 1, 0.0, 0.0);
223     prevFrame.frame->clearValues();
224     ++nextIndex_;
225 }
226
227
228 void
229 AnalysisDataStorage::Impl::notifyPointSet(const AnalysisDataPointSetRef &points)
230 {
231     data_->notifyPointsAdd(points);
232 }
233
234
235 void
236 AnalysisDataStorage::Impl::notifyNextFrames(size_t firstLocation)
237 {
238     if (firstLocation != firstFrameLocation_)
239     {
240         // firstLocation can only be zero here if !storeAll() because
241         // firstFrameLocation_ is always zero for storeAll()
242         int prevIndex =
243             (firstLocation == 0 ? frames_.size() - 1 : firstLocation - 1);
244         if (!frames_[prevIndex].isNotified())
245         {
246             return;
247         }
248     }
249     size_t i = firstLocation;
250     size_t end = endStorageLocation();
251     while (i != end)
252     {
253         Impl::StoredFrame &storedFrame = frames_[i];
254         if (!storedFrame.isFinished())
255         {
256             break;
257         }
258         if (storedFrame.status == StoredFrame::eFinished)
259         {
260             data_->notifyFrameStart(storedFrame.frame->header());
261             data_->notifyPointsAdd(storedFrame.frame->currentPoints());
262             data_->notifyFrameFinish(storedFrame.frame->header());
263             storedFrame.status = StoredFrame::eNotified;
264             if (storedFrame.frame->frameIndex() >= storageLimit_)
265             {
266                 rotateBuffer();
267             }
268         }
269         ++i;
270         if (!storeAll() && i >= frames_.size())
271         {
272             i = 0;
273         }
274     }
275 }
276
277
278 /********************************************************************
279  * AnalysisDataStorage
280  */
281
282 AnalysisDataStorage::AnalysisDataStorage()
283     : impl_(new Impl())
284 {
285 }
286
287
288 AnalysisDataStorage::~AnalysisDataStorage()
289 {
290 }
291
292
293 void
294 AnalysisDataStorage::setMultipoint(bool bMultipoint)
295 {
296     if (bMultipoint && impl_->storageLimit_ > 0)
297     {
298         GMX_THROW(APIError("Storage of multipoint data not supported"));
299     }
300     impl_->bMultipoint_ = bMultipoint;
301 }
302
303
304 void
305 AnalysisDataStorage::setParallelOptions(const AnalysisDataParallelOptions &opt)
306 {
307     impl_->pendingLimit_ = 2 * opt.parallelizationFactor() - 1;
308 }
309
310
311 AnalysisDataFrameRef
312 AnalysisDataStorage::tryGetDataFrame(int index) const
313 {
314     if (impl_->isMultipoint())
315     {
316         return AnalysisDataFrameRef();
317     }
318     int storageIndex = impl_->computeStorageLocation(index);
319     if (storageIndex == -1)
320     {
321         return AnalysisDataFrameRef();
322     }
323     const Impl::StoredFrame &storedFrame = impl_->frames_[storageIndex];
324     if (!storedFrame.isAvailable())
325     {
326         return AnalysisDataFrameRef();
327     }
328     const Impl::FramePointer &frame = storedFrame.frame;
329     return AnalysisDataFrameRef(frame->header(), frame->values_);
330 }
331
332
333 bool
334 AnalysisDataStorage::requestStorage(int nframes)
335 {
336     if (impl_->isMultipoint())
337     {
338         return false;
339     }
340
341     // Handle the case when everything needs to be stored.
342     if (nframes == -1)
343     {
344         impl_->storageLimit_ = std::numeric_limits<int>::max();
345         return true;
346     }
347     // Check whether an earlier call has requested more storage.
348     if (nframes < impl_->storageLimit_)
349     {
350         return true;
351     }
352     impl_->storageLimit_ = nframes;
353     return true;
354 }
355
356
357 void
358 AnalysisDataStorage::startDataStorage(AbstractAnalysisData *data)
359 {
360     // Data needs to be set before calling extendBuffer()
361     impl_->data_ = data;
362     setMultipoint(data->isMultipoint());
363     if (!impl_->storeAll())
364     {
365         impl_->extendBuffer(this, impl_->storageLimit_ + impl_->pendingLimit_ + 1);
366     }
367 }
368
369
370 AnalysisDataStorageFrame &
371 AnalysisDataStorage::startFrame(const AnalysisDataFrameHeader &header)
372 {
373     GMX_ASSERT(header.isValid(), "Invalid header");
374     Impl::StoredFrame *storedFrame;
375     if (impl_->storeAll())
376     {
377         size_t size = header.index() + 1;
378         if (impl_->frames_.size() < size)
379         {
380             impl_->extendBuffer(this, size);
381         }
382         storedFrame = &impl_->frames_[header.index()];
383     }
384     else
385     {
386         int storageIndex = impl_->computeStorageLocation(header.index());
387         if (storageIndex == -1)
388         {
389             GMX_THROW(APIError("Out of bounds frame index"));
390         }
391         storedFrame = &impl_->frames_[storageIndex];
392     }
393     GMX_RELEASE_ASSERT(!storedFrame->isStarted(),
394                        "startFrame() called twice for the same frame");
395     GMX_RELEASE_ASSERT(storedFrame->frame->frameIndex() == header.index(),
396                        "Inconsistent internal frame indexing");
397     storedFrame->status = Impl::StoredFrame::eStarted;
398     storedFrame->frame->header_ = header;
399     if (impl_->isMultipoint())
400     {
401         impl_->data_->notifyFrameStart(header);
402     }
403     return *storedFrame->frame;
404 }
405
406
407 AnalysisDataStorageFrame &
408 AnalysisDataStorage::startFrame(int index, real x, real dx)
409 {
410     return startFrame(AnalysisDataFrameHeader(index, x, dx));
411 }
412
413
414 AnalysisDataStorageFrame &
415 AnalysisDataStorage::currentFrame(int index)
416 {
417     int storageIndex = impl_->computeStorageLocation(index);
418     GMX_RELEASE_ASSERT(storageIndex >= 0, "Out of bounds frame index");
419     Impl::StoredFrame &storedFrame = impl_->frames_[storageIndex];
420     GMX_RELEASE_ASSERT(storedFrame.isStarted(),
421                        "currentFrame() called for frame before startFrame()");
422     GMX_RELEASE_ASSERT(!storedFrame.isFinished(),
423                        "currentFrame() called for frame after finishFrame()");
424     GMX_RELEASE_ASSERT(storedFrame.frame->frameIndex() == index,
425                        "Inconsistent internal frame indexing");
426     return *storedFrame.frame;
427 }
428
429
430 void
431 AnalysisDataStorage::finishFrame(int index)
432 {
433     int storageIndex = impl_->computeStorageLocation(index);
434     GMX_RELEASE_ASSERT(storageIndex >= 0, "Out of bounds frame index");
435     Impl::StoredFrame &storedFrame = impl_->frames_[storageIndex];
436     GMX_RELEASE_ASSERT(storedFrame.isStarted(),
437                        "finishFrame() called for frame before startFrame()");
438     GMX_RELEASE_ASSERT(!storedFrame.isFinished(),
439                        "finishFrame() called twice for the same frame");
440     GMX_RELEASE_ASSERT(storedFrame.frame->frameIndex() == index,
441                        "Inconsistent internal frame indexing");
442     storedFrame.status = Impl::StoredFrame::eFinished;
443     if (impl_->isMultipoint())
444     {
445         // TODO: Check that the last point set has been finished
446         impl_->data_->notifyFrameFinish(storedFrame.frame->header());
447         if (storedFrame.frame->frameIndex() >= impl_->storageLimit_)
448         {
449             impl_->rotateBuffer();
450         }
451     }
452     else
453     {
454         impl_->notifyNextFrames(storageIndex);
455     }
456 }
457
458
459 void
460 AnalysisDataStorage::finishFrame(const AnalysisDataStorageFrame &frame)
461 {
462     finishFrame(frame.frameIndex());
463 }
464
465 } // namespace gmx