Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / analysisdata / datastorage.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013, by the GROMACS development team, led by
5  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6  * others, as listed in the AUTHORS file in the top-level source
7  * directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS 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 GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements classes in datastorage.h and paralleloptions.h.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_analysisdata
41  */
42 #include "datastorage.h"
43
44 #include <algorithm>
45 #include <iterator>
46 #include <limits>
47 #include <vector>
48
49 #include "gromacs/analysisdata/abstractdata.h"
50 #include "gromacs/analysisdata/dataframe.h"
51 #include "gromacs/analysisdata/datamodulemanager.h"
52 #include "gromacs/analysisdata/paralleloptions.h"
53 #include "gromacs/utility/exceptions.h"
54 #include "gromacs/utility/gmxassert.h"
55 #include "gromacs/utility/uniqueptr.h"
56
57 namespace gmx
58 {
59
60 /********************************************************************
61  * AnalysisDataParallelOptions
62  */
63
64 AnalysisDataParallelOptions::AnalysisDataParallelOptions()
65     : parallelizationFactor_(1)
66 {
67 }
68
69
70 AnalysisDataParallelOptions::AnalysisDataParallelOptions(int parallelizationFactor)
71     : parallelizationFactor_(parallelizationFactor)
72 {
73     GMX_RELEASE_ASSERT(parallelizationFactor >= 1,
74                        "Invalid parallelization factor");
75 }
76
77
78 /********************************************************************
79  * AnalysisDataStorageImpl declaration
80  */
81
82 namespace internal
83 {
84
85 //! Smart pointer type for managing a storage frame builder.
86 typedef gmx_unique_ptr<AnalysisDataStorageFrame>::type
87     AnalysisDataFrameBuilderPointer;
88
89 /*! \internal \brief
90  * Private implementation class for AnalysisDataStorage.
91  *
92  * \ingroup module_analysisdata
93  */
94 class AnalysisDataStorageImpl
95 {
96     public:
97         //! Smart pointer type for managing a stored frame.
98         typedef gmx_unique_ptr<AnalysisDataStorageFrameData>::type FramePointer;
99
100         //! Shorthand for a list of data frames that are currently stored.
101         typedef std::vector<FramePointer> FrameList;
102         //! Shorthand for a list of currently unused storage frame builders.
103         typedef std::vector<AnalysisDataFrameBuilderPointer> FrameBuilderList;
104
105         AnalysisDataStorageImpl();
106
107         //! Returns whether the storage is set to use multipoint data.
108         bool isMultipoint() const;
109         /*! \brief
110          * Whether storage of all frames has been requested.
111          *
112          * Storage of all frames also works as expected if \a storageLimit_ is
113          * used in comparisons directly, but this method should be used to
114          * check how to manage \a frames_.
115          */
116         bool storeAll() const
117         {
118             return storageLimit_ == std::numeric_limits<int>::max();
119         }
120         //! Returns the index of the oldest frame that may be currently stored.
121         int firstStoredIndex() const;
122         //! Returns the index of the first frame that is not fully notified.
123         int firstUnnotifiedIndex() const { return firstUnnotifiedIndex_; }
124         /*! \brief
125          * Computes index into \a frames_ for accessing frame \p index.
126          *
127          * \param[in]  index  Zero-based frame index.
128          * \retval  -1 if \p index is not available in \a frames_.
129          *
130          * Does not throw.
131          */
132         int computeStorageLocation(int index) const;
133
134         /*! \brief
135          * Computes an index into \a frames_ that is one past the last frame
136          * stored.
137          *
138          * Does not throw.
139          */
140         size_t endStorageLocation() const;
141
142         /*! \brief
143          * Extends \a frames_ to a new size.
144          *
145          * \throws std::bad_alloc if out of memory.
146          */
147         void extendBuffer(size_t newSize);
148         /*! \brief
149          * Remove oldest frame from the storage to make space for a new one.
150          *
151          * Increments \a firstFrameLocation_ and reinitializes the frame that
152          * was made unavailable by this operation.
153          *
154          * Does not throw.
155          *
156          * \see frames_
157          */
158         void rotateBuffer();
159
160         /*! \brief
161          * Returns a frame builder object for use with a new frame.
162          *
163          * \throws std::bad_alloc if out of memory.
164          */
165         AnalysisDataFrameBuilderPointer getFrameBuilder();
166
167         /*! \brief
168          * Returns whether notifications should be immediately fired.
169          *
170          * This is used to optimize multipoint handling for non-parallel cases,
171          * where it is not necessary to store even a single frame.
172          *
173          * Does not throw.
174          */
175         bool shouldNotifyImmediately() const
176         {
177             return isMultipoint() && storageLimit_ == 0 && pendingLimit_ == 1;
178         }
179         /*! \brief
180          * Returns whether data needs to be stored at all.
181          *
182          * This is used to optimize multipoint handling for parallel cases
183          * (where shouldNotifyImmediately() returns false),
184          * where it is not necessary to store even a single frame.
185          *
186          * \todo
187          * This could be extended to non-multipoint data as well.
188          *
189          * Does not throw.
190          */
191         bool needStorage() const
192         {
193             return storageLimit_ > 0 || (pendingLimit_ > 1 && modules_->hasSerialModules());
194         }
195         /*! \brief
196          * Calls notification methods for new frames.
197          *
198          * \param[in] firstLocation  First frame to consider.
199          * \throws    unspecified  Any exception thrown by frame notification
200          *      methods in AbstractAnalysisData.
201          *
202          * Notifies \a data_ of new frames (from \p firstLocation and after
203          * that) if all previous frames have already been notified.
204          * Also rotates the \a frames_ buffer as necessary.
205          */
206         void notifyNextFrames(size_t firstLocation);
207         //! Implementation for AnalysisDataStorage::finishFrame().
208         void finishFrame(int index);
209
210
211         //! Parent data object to access data dimensionality etc.
212         const AbstractAnalysisData *data_;
213         //! Manager to use for notification calls.
214         AnalysisDataModuleManager  *modules_;
215         /*! \brief
216          * Number of past frames that need to be stored.
217          *
218          * Always non-negative.  If storage of all frames has been requested,
219          * this is set to a large number.
220          */
221         int                     storageLimit_;
222         /*! \brief
223          * Number of future frames that may need to be started.
224          *
225          * Should always be at least one.
226          *
227          * \todo
228          * Get rid of this alltogether, as it is no longer used much.
229          *
230          * \see AnalysisDataStorage::startFrame()
231          */
232         int                     pendingLimit_;
233         /*! \brief
234          * Data frames that are currently stored.
235          *
236          * If storage of all frames has been requested, this is simply a vector
237          * of frames up to the latest frame that has been started.
238          * In this case, \a firstFrameLocation_ is always zero.
239          *
240          * If storage of all frames is not requested, this is a ring buffer of
241          * frames of size \c n=storageLimit_+pendingLimit_+1.  If a frame with
242          * index \c index is currently stored, its location is
243          * \c index%frames_.size().
244          * When at most \a storageLimit_ first frames have been finished,
245          * this contains storage for the first \c n-1 frames.
246          * When more than \a storageLimit_ first frames have been finished,
247          * the oldest stored frame is stored in the location
248          * \a firstFrameLocation_, and \a storageLimit_ frames starting from
249          * this location are the last finished frames.  \a pendingLimit_ frames
250          * follow, and some of these may be in progress or finished.
251          * There is always one unused frame in the buffer, which is initialized
252          * such that when \a firstFrameLocation_ is incremented, it becomes
253          * valid.  This makes it easier to rotate the buffer in concurrent
254          * access scenarions (which are not yet otherwise implemented).
255          */
256         FrameList               frames_;
257         //! Location of oldest frame in \a frames_.
258         size_t                  firstFrameLocation_;
259         //! Index of the first frame that is not fully notified.
260         int                     firstUnnotifiedIndex_;
261         /*! \brief
262          * Currently unused frame builders.
263          *
264          * The builders are cached to avoid repeatedly allocating memory for
265          * them.  Typically, there are as many builders as there are concurrent
266          * users of the storage object.  Whenever a frame is started, a builder
267          * is pulled from this pool by getFrameBuilder() (a new one is created
268          * if none are available), and assigned for that frame.  When that
269          * frame is finished, the builder is returned to this pool.
270          */
271         FrameBuilderList        builders_;
272         /*! \brief
273          * Index of next frame that will be added to \a frames_.
274          *
275          * If all frames are not stored, this will be the index of the unused
276          * frame (see \a frames_).
277          */
278         int                     nextIndex_;
279 };
280
281 /********************************************************************
282  * AnalysisDataStorageFrameImpl declaration
283  */
284
285 /*! \internal \brief
286  * Internal representation for a single stored frame.
287  *
288  * It is implemented such that the frame header is always valid, i.e.,
289  * header().isValid() returns always true.
290  *
291  * Methods in this class do not throw unless otherwise indicated.
292  *
293  * \ingroup module_analysisdata
294  */
295 class AnalysisDataStorageFrameData
296 {
297     public:
298         //! Shorthand for a iterator into storage value containers.
299         typedef std::vector<AnalysisDataValue>::const_iterator ValueIterator;
300
301         //! Indicates what operations have been performed on a frame.
302         enum Status
303         {
304             eMissing,  //!< Frame has not yet been started.
305             eStarted,  //!< startFrame() has been called.
306             eFinished, //!< finishFrame() has been called.
307             eNotified  //!< Appropriate notifications have been sent.
308         };
309
310         /*! \brief
311          * Create a new storage frame.
312          *
313          * \param     storageImpl  Storage object this frame belongs to.
314          * \param[in] index        Zero-based index for the frame.
315          */
316         AnalysisDataStorageFrameData(AnalysisDataStorageImpl *storageImpl,
317                                      int                      index);
318
319         //! Whether the frame has been started with startFrame().
320         bool isStarted() const { return status_ >= eStarted; }
321         //! Whether the frame has been finished with finishFrame().
322         bool isFinished() const { return status_ >= eFinished; }
323         //! Whether all notifications have been sent.
324         bool isNotified() const { return status_ >= eNotified; }
325         //! Whether the frame is ready to be available outside the storage.
326         bool isAvailable() const { return status_ >= eFinished; }
327
328         //! Marks the frame as notified.
329         void markNotified() { status_ = eNotified; }
330
331         //! Returns the storage implementation object.
332         AnalysisDataStorageImpl &storageImpl() const { return storageImpl_; }
333         //! Returns the underlying data object (for data dimensionalities etc.).
334         const AbstractAnalysisData &baseData() const { return *storageImpl().data_; }
335
336         //! Returns header for the frame.
337         const AnalysisDataFrameHeader &header() const { return header_; }
338         //! Returns zero-based index of the frame.
339         int frameIndex() const { return header().index(); }
340         //! Returns the number of point sets for the frame.
341         int pointSetCount() const { return pointSets_.size(); }
342
343         //! Clears the frame for reusing as a new frame.
344         void clearFrame(int newIndex);
345         /*! \brief
346          * Initializes the frame during AnalysisDataStorage::startFrame().
347          *
348          * \param[in] header  Header to use for the new frame.
349          * \param[in] builder Builder object to use.
350          */
351         void startFrame(const AnalysisDataFrameHeader   &header,
352                         AnalysisDataFrameBuilderPointer  builder);
353         //! Returns the builder for this frame.
354         AnalysisDataStorageFrame &builder() const
355         {
356             GMX_ASSERT(builder_, "Accessing builder for not-in-progress frame");
357             return *builder_;
358         }
359         /*! \brief
360          * Adds a new point set to this frame.
361          */
362         void addPointSet(int dataSetIndex, int firstColumn,
363                          ValueIterator begin, ValueIterator end);
364         /*! \brief
365          * Finalizes the frame during AnalysisDataStorage::finishFrame().
366          *
367          * \returns The builder object used by the frame, for reusing it for
368          *      other frames.
369          */
370         AnalysisDataFrameBuilderPointer finishFrame(bool bMultipoint);
371
372         //! Returns frame reference to this frame.
373         AnalysisDataFrameRef frameReference() const
374         {
375             return AnalysisDataFrameRef(header_, values_, pointSets_);
376         }
377         //! Returns point set reference to a given point set.
378         AnalysisDataPointSetRef pointSet(int index) const;
379
380     private:
381         //! Storage object that contains this frame.
382         AnalysisDataStorageImpl                &storageImpl_;
383         //! Header for the frame.
384         AnalysisDataFrameHeader                 header_;
385         //! Values for the frame.
386         std::vector<AnalysisDataValue>          values_;
387         //! Information about each point set in the frame.
388         std::vector<AnalysisDataPointSetInfo>   pointSets_;
389         /*! \brief
390          * Builder object for the frame.
391          *
392          * Non-NULL when the frame is in progress, i.e., has been started but
393          * not yet finished.
394          */
395         AnalysisDataFrameBuilderPointer         builder_;
396         //! In what state the frame currently is.
397         Status                                  status_;
398
399         GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisDataStorageFrameData);
400 };
401
402 /********************************************************************
403  * AnalysisDataStorageImpl implementation
404  */
405
406 AnalysisDataStorageImpl::AnalysisDataStorageImpl()
407     : data_(NULL), modules_(NULL),
408       storageLimit_(0), pendingLimit_(1),
409       firstFrameLocation_(0), firstUnnotifiedIndex_(0), nextIndex_(0)
410 {
411 }
412
413
414 bool
415 AnalysisDataStorageImpl::isMultipoint() const
416 {
417     GMX_ASSERT(data_ != NULL, "isMultipoint() called too early");
418     return data_->isMultipoint();
419 }
420
421
422 int
423 AnalysisDataStorageImpl::firstStoredIndex() const
424 {
425     return frames_[firstFrameLocation_]->frameIndex();
426 }
427
428
429 int
430 AnalysisDataStorageImpl::computeStorageLocation(int index) const
431 {
432     if (index < firstStoredIndex() || index >= nextIndex_)
433     {
434         return -1;
435     }
436     return index % frames_.size();
437 }
438
439
440 size_t
441 AnalysisDataStorageImpl::endStorageLocation() const
442 {
443     if (storeAll())
444     {
445         return frames_.size();
446     }
447     if (frames_[0]->frameIndex() == 0 || firstFrameLocation_ == 0)
448     {
449         return frames_.size() - 1;
450     }
451     return firstFrameLocation_ - 1;
452 }
453
454
455 void
456 AnalysisDataStorageImpl::extendBuffer(size_t newSize)
457 {
458     frames_.reserve(newSize);
459     while (frames_.size() < newSize)
460     {
461         frames_.push_back(FramePointer(new AnalysisDataStorageFrameData(this, nextIndex_)));
462         ++nextIndex_;
463     }
464     // The unused frame should not be included in the count.
465     if (!storeAll())
466     {
467         --nextIndex_;
468     }
469 }
470
471
472 void
473 AnalysisDataStorageImpl::rotateBuffer()
474 {
475     GMX_ASSERT(!storeAll(),
476                "No need to rotate internal buffer if everything is stored");
477     size_t prevFirst = firstFrameLocation_;
478     size_t nextFirst = prevFirst + 1;
479     if (nextFirst == frames_.size())
480     {
481         nextFirst = 0;
482     }
483     firstFrameLocation_ = nextFirst;
484     frames_[prevFirst]->clearFrame(nextIndex_ + 1);
485     ++nextIndex_;
486 }
487
488
489 AnalysisDataFrameBuilderPointer
490 AnalysisDataStorageImpl::getFrameBuilder()
491 {
492     if (builders_.empty())
493     {
494         return AnalysisDataFrameBuilderPointer(new AnalysisDataStorageFrame(*data_));
495     }
496     AnalysisDataFrameBuilderPointer builder(move(builders_.back()));
497     builders_.pop_back();
498     return move(builder);
499 }
500
501
502 void
503 AnalysisDataStorageImpl::notifyNextFrames(size_t firstLocation)
504 {
505     if (frames_[firstLocation]->frameIndex() != firstUnnotifiedIndex_)
506     {
507         return;
508     }
509     size_t i   = firstLocation;
510     size_t end = endStorageLocation();
511     while (i != end)
512     {
513         AnalysisDataStorageFrameData &storedFrame = *frames_[i];
514         if (!storedFrame.isFinished())
515         {
516             break;
517         }
518         if (!storedFrame.isNotified())
519         {
520             // Increment before the notifications to make the frame available
521             // in the module callbacks.
522             ++firstUnnotifiedIndex_;
523             modules_->notifyFrameStart(storedFrame.header());
524             for (int j = 0; j < storedFrame.pointSetCount(); ++j)
525             {
526                 modules_->notifyPointsAdd(storedFrame.pointSet(j));
527             }
528             modules_->notifyFrameFinish(storedFrame.header());
529             storedFrame.markNotified();
530             if (storedFrame.frameIndex() >= storageLimit_)
531             {
532                 rotateBuffer();
533             }
534         }
535         ++i;
536         if (!storeAll() && i >= frames_.size())
537         {
538             i = 0;
539         }
540     }
541 }
542
543
544 void
545 AnalysisDataStorageImpl::finishFrame(int index)
546 {
547     const int storageIndex = computeStorageLocation(index);
548     GMX_RELEASE_ASSERT(storageIndex >= 0, "Out of bounds frame index");
549
550     AnalysisDataStorageFrameData &storedFrame = *frames_[storageIndex];
551     GMX_RELEASE_ASSERT(storedFrame.isStarted(),
552                        "finishFrame() called for frame before startFrame()");
553     GMX_RELEASE_ASSERT(!storedFrame.isFinished(),
554                        "finishFrame() called twice for the same frame");
555     GMX_RELEASE_ASSERT(storedFrame.frameIndex() == index,
556                        "Inconsistent internal frame indexing");
557     builders_.push_back(storedFrame.finishFrame(isMultipoint()));
558     modules_->notifyParallelFrameFinish(storedFrame.header());
559     if (shouldNotifyImmediately())
560     {
561         ++firstUnnotifiedIndex_;
562         modules_->notifyFrameFinish(storedFrame.header());
563         if (storedFrame.frameIndex() >= storageLimit_)
564         {
565             rotateBuffer();
566         }
567     }
568     else
569     {
570         notifyNextFrames(storageIndex);
571     }
572 }
573
574
575 /********************************************************************
576  * AnalysisDataStorageFrame implementation
577  */
578
579 AnalysisDataStorageFrameData::AnalysisDataStorageFrameData(
580         AnalysisDataStorageImpl *storageImpl,
581         int                      index)
582     : storageImpl_(*storageImpl), header_(index, 0.0, 0.0), status_(eMissing)
583 {
584     GMX_RELEASE_ASSERT(storageImpl->data_ != NULL,
585                        "Storage frame constructed before data started");
586     // With non-multipoint data, the point set structure is static,
587     // so initialize it only once here.
588     if (!baseData().isMultipoint())
589     {
590         int offset = 0;
591         for (int i = 0; i < baseData().dataSetCount(); ++i)
592         {
593             int columnCount = baseData().columnCount(i);
594             pointSets_.push_back(
595                     AnalysisDataPointSetInfo(offset, columnCount, i, 0));
596             offset += columnCount;
597         }
598     }
599 }
600
601
602 void
603 AnalysisDataStorageFrameData::clearFrame(int newIndex)
604 {
605     GMX_RELEASE_ASSERT(!builder_, "Should not clear an in-progress frame");
606     status_ = eMissing;
607     header_ = AnalysisDataFrameHeader(newIndex, 0.0, 0.0);
608     values_.clear();
609     if (baseData().isMultipoint())
610     {
611         pointSets_.clear();
612     }
613 }
614
615
616 void
617 AnalysisDataStorageFrameData::startFrame(
618         const AnalysisDataFrameHeader   &header,
619         AnalysisDataFrameBuilderPointer  builder)
620 {
621     status_         = eStarted;
622     header_         = header;
623     builder_        = move(builder);
624     builder_->data_ = this;
625     builder_->selectDataSet(0);
626 }
627
628
629 void
630 AnalysisDataStorageFrameData::addPointSet(int dataSetIndex, int firstColumn,
631                                           ValueIterator begin, ValueIterator end)
632 {
633     const int                valueCount = end - begin;
634     AnalysisDataPointSetInfo pointSetInfo(0, valueCount,
635                                           dataSetIndex, firstColumn);
636     AnalysisDataPointSetRef  pointSet(header(), pointSetInfo,
637                                       AnalysisDataValuesRef(begin, end));
638     storageImpl().modules_->notifyParallelPointsAdd(pointSet);
639     if (storageImpl().shouldNotifyImmediately())
640     {
641         storageImpl().modules_->notifyPointsAdd(pointSet);
642     }
643     else if (storageImpl().needStorage())
644     {
645         pointSets_.push_back(
646                 AnalysisDataPointSetInfo(values_.size(), valueCount,
647                                          dataSetIndex, firstColumn));
648         std::copy(begin, end, std::back_inserter(values_));
649     }
650 }
651
652
653 AnalysisDataFrameBuilderPointer
654 AnalysisDataStorageFrameData::finishFrame(bool bMultipoint)
655 {
656     status_ = eFinished;
657     if (!bMultipoint)
658     {
659         GMX_RELEASE_ASSERT(static_cast<int>(pointSets_.size()) == baseData().dataSetCount(),
660                            "Point sets created for non-multipoint data");
661         values_ = builder_->values_;
662         builder_->clearValues();
663         for (int i = 0; i < pointSetCount(); ++i)
664         {
665             storageImpl().modules_->notifyParallelPointsAdd(pointSet(i));
666         }
667     }
668     else
669     {
670         GMX_RELEASE_ASSERT(!builder_->bPointSetInProgress_,
671                            "Unfinished point set");
672     }
673     AnalysisDataFrameBuilderPointer builder(move(builder_));
674     builder_.reset();
675     return move(builder);
676 }
677
678
679 AnalysisDataPointSetRef
680 AnalysisDataStorageFrameData::pointSet(int index) const
681 {
682     GMX_ASSERT(index >= 0 && index < pointSetCount(),
683                "Invalid point set index");
684     return AnalysisDataPointSetRef(
685             header_, pointSets_[index],
686             AnalysisDataValuesRef(values_.begin(), values_.end()));
687 }
688
689 }   // namespace internal
690
691 /********************************************************************
692  * AnalysisDataStorageFrame
693  */
694
695 AnalysisDataStorageFrame::AnalysisDataStorageFrame(
696         const AbstractAnalysisData &data)
697     : data_(NULL), currentDataSet_(0), currentOffset_(0),
698       columnCount_(data.columnCount(0)), bPointSetInProgress_(false)
699 {
700     int totalColumnCount = 0;
701     for (int i = 0; i < data.dataSetCount(); ++i)
702     {
703         totalColumnCount += data.columnCount(i);
704     }
705     values_.resize(totalColumnCount);
706 }
707
708
709 AnalysisDataStorageFrame::~AnalysisDataStorageFrame()
710 {
711 }
712
713
714 void
715 AnalysisDataStorageFrame::clearValues()
716 {
717     if (bPointSetInProgress_)
718     {
719         std::vector<AnalysisDataValue>::iterator i;
720         for (i = values_.begin(); i != values_.end(); ++i)
721         {
722             i->clear();
723         }
724     }
725     bPointSetInProgress_ = false;
726 }
727
728
729 void
730 AnalysisDataStorageFrame::selectDataSet(int index)
731 {
732     GMX_RELEASE_ASSERT(data_ != NULL, "Invalid frame accessed");
733     const AbstractAnalysisData &baseData = data_->baseData();
734     GMX_RELEASE_ASSERT(index >= 0 && index < baseData.dataSetCount(),
735                        "Out of range data set index");
736     GMX_RELEASE_ASSERT(!baseData.isMultipoint() || !bPointSetInProgress_,
737                        "Point sets in multipoint data cannot span data sets");
738     currentDataSet_ = index;
739     currentOffset_  = 0;
740     // TODO: Consider precalculating.
741     for (int i = 0; i < index; ++i)
742     {
743         currentOffset_ += baseData.columnCount(i);
744     }
745     columnCount_    = baseData.columnCount(index);
746 }
747
748
749 void
750 AnalysisDataStorageFrame::finishPointSet()
751 {
752     GMX_RELEASE_ASSERT(data_ != NULL, "Invalid frame accessed");
753     GMX_RELEASE_ASSERT(data_->baseData().isMultipoint(),
754                        "Should not be called for non-multipoint data");
755     if (bPointSetInProgress_)
756     {
757         std::vector<AnalysisDataValue>::const_iterator begin
758             = values_.begin() + currentOffset_;
759         std::vector<AnalysisDataValue>::const_iterator end
760             = begin + columnCount_;
761         int firstColumn = 0;
762         while (begin != end && !begin->isSet())
763         {
764             ++begin;
765             ++firstColumn;
766         }
767         while (end != begin && !(end-1)->isSet())
768         {
769             --end;
770         }
771         if (begin == end)
772         {
773             firstColumn = 0;
774         }
775         data_->addPointSet(currentDataSet_, firstColumn, begin, end);
776     }
777     clearValues();
778 }
779
780
781 void
782 AnalysisDataStorageFrame::finishFrame()
783 {
784     GMX_RELEASE_ASSERT(data_ != NULL, "Invalid frame accessed");
785     data_->storageImpl().finishFrame(data_->frameIndex());
786 }
787
788
789 /********************************************************************
790  * AnalysisDataStorage
791  */
792
793 AnalysisDataStorage::AnalysisDataStorage()
794     : impl_(new Impl())
795 {
796 }
797
798
799 AnalysisDataStorage::~AnalysisDataStorage()
800 {
801 }
802
803
804 int
805 AnalysisDataStorage::frameCount() const
806 {
807     return impl_->firstUnnotifiedIndex();
808 }
809
810
811 AnalysisDataFrameRef
812 AnalysisDataStorage::tryGetDataFrame(int index) const
813 {
814     int storageIndex = impl_->computeStorageLocation(index);
815     if (storageIndex == -1)
816     {
817         return AnalysisDataFrameRef();
818     }
819     const internal::AnalysisDataStorageFrameData &storedFrame
820         = *impl_->frames_[storageIndex];
821     if (!storedFrame.isAvailable())
822     {
823         return AnalysisDataFrameRef();
824     }
825     return storedFrame.frameReference();
826 }
827
828
829 bool
830 AnalysisDataStorage::requestStorage(int nframes)
831 {
832     // Handle the case when everything needs to be stored.
833     if (nframes == -1)
834     {
835         impl_->storageLimit_ = std::numeric_limits<int>::max();
836         return true;
837     }
838     // Check whether an earlier call has requested more storage.
839     if (nframes < impl_->storageLimit_)
840     {
841         return true;
842     }
843     impl_->storageLimit_ = nframes;
844     return true;
845 }
846
847
848 void
849 AnalysisDataStorage::startDataStorage(AbstractAnalysisData      *data,
850                                       AnalysisDataModuleManager *modules)
851 {
852     modules->notifyDataStart(data);
853     // Data needs to be set before calling extendBuffer()
854     impl_->data_    = data;
855     impl_->modules_ = modules;
856     if (!impl_->storeAll())
857     {
858         // 2 = pending limit (1) + 1
859         impl_->extendBuffer(impl_->storageLimit_ + 2);
860     }
861 }
862
863
864 void
865 AnalysisDataStorage::startParallelDataStorage(
866         AbstractAnalysisData              *data,
867         AnalysisDataModuleManager         *modules,
868         const AnalysisDataParallelOptions &options)
869 {
870     const int pendingLimit = 2 * options.parallelizationFactor() - 1;
871     impl_->pendingLimit_   = pendingLimit;
872     modules->notifyParallelDataStart(data, options);
873     // Data needs to be set before calling extendBuffer()
874     impl_->data_    = data;
875     impl_->modules_ = modules;
876     if (!impl_->storeAll())
877     {
878         impl_->extendBuffer(impl_->storageLimit_ + pendingLimit + 1);
879     }
880 }
881
882
883 AnalysisDataStorageFrame &
884 AnalysisDataStorage::startFrame(const AnalysisDataFrameHeader &header)
885 {
886     GMX_ASSERT(header.isValid(), "Invalid header");
887     internal::AnalysisDataStorageFrameData *storedFrame;
888     if (impl_->storeAll())
889     {
890         size_t size = header.index() + 1;
891         if (impl_->frames_.size() < size)
892         {
893             impl_->extendBuffer(size);
894         }
895         storedFrame = impl_->frames_[header.index()].get();
896     }
897     else
898     {
899         int storageIndex = impl_->computeStorageLocation(header.index());
900         if (storageIndex == -1)
901         {
902             GMX_THROW(APIError("Out of bounds frame index"));
903         }
904         storedFrame = impl_->frames_[storageIndex].get();
905     }
906     GMX_RELEASE_ASSERT(!storedFrame->isStarted(),
907                        "startFrame() called twice for the same frame");
908     GMX_RELEASE_ASSERT(storedFrame->frameIndex() == header.index(),
909                        "Inconsistent internal frame indexing");
910     storedFrame->startFrame(header, impl_->getFrameBuilder());
911     impl_->modules_->notifyParallelFrameStart(header);
912     if (impl_->shouldNotifyImmediately())
913     {
914         impl_->modules_->notifyFrameStart(header);
915     }
916     return storedFrame->builder();
917 }
918
919
920 AnalysisDataStorageFrame &
921 AnalysisDataStorage::startFrame(int index, real x, real dx)
922 {
923     return startFrame(AnalysisDataFrameHeader(index, x, dx));
924 }
925
926
927 AnalysisDataStorageFrame &
928 AnalysisDataStorage::currentFrame(int index)
929 {
930     const int storageIndex = impl_->computeStorageLocation(index);
931     GMX_RELEASE_ASSERT(storageIndex >= 0, "Out of bounds frame index");
932
933     internal::AnalysisDataStorageFrameData &storedFrame = *impl_->frames_[storageIndex];
934     GMX_RELEASE_ASSERT(storedFrame.isStarted(),
935                        "currentFrame() called for frame before startFrame()");
936     GMX_RELEASE_ASSERT(!storedFrame.isFinished(),
937                        "currentFrame() called for frame after finishFrame()");
938     GMX_RELEASE_ASSERT(storedFrame.frameIndex() == index,
939                        "Inconsistent internal frame indexing");
940     return storedFrame.builder();
941 }
942
943
944 void
945 AnalysisDataStorage::finishFrame(int index)
946 {
947     impl_->finishFrame(index);
948 }
949
950 void
951 AnalysisDataStorage::finishDataStorage()
952 {
953     // TODO: Check that all frames have been finished etc.
954     impl_->builders_.clear();
955     impl_->modules_->notifyDataFinish();
956 }
957
958 } // namespace gmx