Use smart pointers for internal memory management.
[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/fatalerror/exceptions.h"
46 #include "gromacs/fatalerror/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         AnalysisDataStorageFrame *frame =
197             new AnalysisDataStorageFrame(storage, columnCount(), nextIndex_);
198         frames_.push_back(StoredFrame(frame));
199         ++nextIndex_;
200     }
201     // The unused frame should not be included in the count.
202     if (!storeAll())
203     {
204         --nextIndex_;
205     }
206 }
207
208
209 void
210 AnalysisDataStorage::Impl::rotateBuffer()
211 {
212     GMX_ASSERT(!storeAll(),
213                "No need to rotate internal buffer if everything is stored");
214     size_t prevFirst = firstFrameLocation_;
215     size_t nextFirst = prevFirst + 1;
216     if (nextFirst == frames_.size())
217     {
218         nextFirst = 0;
219     }
220     firstFrameLocation_ = nextFirst;
221     StoredFrame &prevFrame = frames_[prevFirst];
222     prevFrame.status = StoredFrame::eMissing;
223     prevFrame.frame->header_ = AnalysisDataFrameHeader(nextIndex_ + 1, 0.0, 0.0);
224     prevFrame.frame->clearValues();
225     ++nextIndex_;
226 }
227
228
229 void
230 AnalysisDataStorage::Impl::notifyPointSet(const AnalysisDataPointSetRef &points)
231 {
232     data_->notifyPointsAdd(points);
233 }
234
235
236 void
237 AnalysisDataStorage::Impl::notifyNextFrames(size_t firstLocation)
238 {
239     if (firstLocation != firstFrameLocation_)
240     {
241         // firstLocation can only be zero here if !storeAll() because
242         // firstFrameLocation_ is always zero for storeAll()
243         int prevIndex =
244             (firstLocation == 0 ? frames_.size() - 1 : firstLocation - 1);
245         if (!frames_[prevIndex].isNotified())
246         {
247             return;
248         }
249     }
250     size_t i = firstLocation;
251     size_t end = endStorageLocation();
252     while (i != end)
253     {
254         Impl::StoredFrame &storedFrame = frames_[i];
255         if (!storedFrame.isFinished())
256         {
257             break;
258         }
259         if (storedFrame.status == StoredFrame::eFinished)
260         {
261             data_->notifyFrameStart(storedFrame.frame->header());
262             data_->notifyPointsAdd(storedFrame.frame->currentPoints());
263             data_->notifyFrameFinish(storedFrame.frame->header());
264             storedFrame.status = StoredFrame::eNotified;
265             if (storedFrame.frame->frameIndex() >= storageLimit_)
266             {
267                 rotateBuffer();
268             }
269         }
270         ++i;
271         if (!storeAll() && i >= frames_.size())
272         {
273             i = 0;
274         }
275     }
276 }
277
278
279 /********************************************************************
280  * AnalysisDataStorage
281  */
282
283 AnalysisDataStorage::AnalysisDataStorage()
284     : impl_(new Impl())
285 {
286 }
287
288
289 AnalysisDataStorage::~AnalysisDataStorage()
290 {
291 }
292
293
294 void
295 AnalysisDataStorage::setMultipoint(bool bMultipoint)
296 {
297     if (bMultipoint && impl_->storageLimit_ > 0)
298     {
299         GMX_THROW(APIError("Storage of multipoint data not supported"));
300     }
301     impl_->bMultipoint_ = bMultipoint;
302 }
303
304
305 void
306 AnalysisDataStorage::setParallelOptions(const AnalysisDataParallelOptions &opt)
307 {
308     impl_->pendingLimit_ = 2 * opt.parallelizationFactor() - 1;
309 }
310
311
312 AnalysisDataFrameRef
313 AnalysisDataStorage::tryGetDataFrame(int index) const
314 {
315     if (impl_->isMultipoint())
316     {
317         return AnalysisDataFrameRef();
318     }
319     int storageIndex = impl_->computeStorageLocation(index);
320     if (storageIndex == -1)
321     {
322         return AnalysisDataFrameRef();
323     }
324     const Impl::StoredFrame &storedFrame = impl_->frames_[storageIndex];
325     if (!storedFrame.isAvailable())
326     {
327         return AnalysisDataFrameRef();
328     }
329     const Impl::FramePointer &frame = storedFrame.frame;
330     return AnalysisDataFrameRef(frame->header(), frame->values_);
331 }
332
333
334 bool
335 AnalysisDataStorage::requestStorage(int nframes)
336 {
337     if (impl_->isMultipoint())
338     {
339         return false;
340     }
341
342     // Handle the case when everything needs to be stored.
343     if (nframes == -1)
344     {
345         impl_->storageLimit_ = std::numeric_limits<int>::max();
346         return true;
347     }
348     // Check whether an earlier call has requested more storage.
349     if (nframes < impl_->storageLimit_)
350     {
351         return true;
352     }
353     impl_->storageLimit_ = nframes;
354     return true;
355 }
356
357
358 void
359 AnalysisDataStorage::startDataStorage(AbstractAnalysisData *data)
360 {
361     // Data needs to be set before calling extendBuffer()
362     impl_->data_ = data;
363     setMultipoint(data->isMultipoint());
364     if (!impl_->storeAll())
365     {
366         impl_->extendBuffer(this, impl_->storageLimit_ + impl_->pendingLimit_ + 1);
367     }
368 }
369
370
371 AnalysisDataStorageFrame &
372 AnalysisDataStorage::startFrame(const AnalysisDataFrameHeader &header)
373 {
374     GMX_ASSERT(header.isValid(), "Invalid header");
375     Impl::StoredFrame *storedFrame;
376     if (impl_->storeAll())
377     {
378         size_t size = header.index() + 1;
379         if (impl_->frames_.size() < size)
380         {
381             impl_->extendBuffer(this, size);
382         }
383         storedFrame = &impl_->frames_[header.index()];
384     }
385     else
386     {
387         int storageIndex = impl_->computeStorageLocation(header.index());
388         if (storageIndex == -1)
389         {
390             GMX_THROW(APIError("Out of bounds frame index"));
391         }
392         storedFrame = &impl_->frames_[storageIndex];
393     }
394     GMX_RELEASE_ASSERT(!storedFrame->isStarted(),
395                        "startFrame() called twice for the same frame");
396     GMX_RELEASE_ASSERT(storedFrame->frame->frameIndex() == header.index(),
397                        "Inconsistent internal frame indexing");
398     storedFrame->status = Impl::StoredFrame::eStarted;
399     storedFrame->frame->header_ = header;
400     if (impl_->isMultipoint())
401     {
402         impl_->data_->notifyFrameStart(header);
403     }
404     return *storedFrame->frame;
405 }
406
407
408 AnalysisDataStorageFrame &
409 AnalysisDataStorage::startFrame(int index, real x, real dx)
410 {
411     return startFrame(AnalysisDataFrameHeader(index, x, dx));
412 }
413
414
415 AnalysisDataStorageFrame &
416 AnalysisDataStorage::currentFrame(int index)
417 {
418     int storageIndex = impl_->computeStorageLocation(index);
419     GMX_RELEASE_ASSERT(storageIndex >= 0, "Out of bounds frame index");
420     Impl::StoredFrame &storedFrame = impl_->frames_[storageIndex];
421     GMX_RELEASE_ASSERT(storedFrame.isStarted(),
422                        "currentFrame() called for frame before startFrame()");
423     GMX_RELEASE_ASSERT(!storedFrame.isFinished(),
424                        "currentFrame() called for frame after finishFrame()");
425     GMX_RELEASE_ASSERT(storedFrame.frame->frameIndex() == index,
426                        "Inconsistent internal frame indexing");
427     return *storedFrame.frame;
428 }
429
430
431 void
432 AnalysisDataStorage::finishFrame(int index)
433 {
434     int storageIndex = impl_->computeStorageLocation(index);
435     GMX_RELEASE_ASSERT(storageIndex >= 0, "Out of bounds frame index");
436     Impl::StoredFrame &storedFrame = impl_->frames_[storageIndex];
437     GMX_RELEASE_ASSERT(storedFrame.isStarted(),
438                        "finishFrame() called for frame before startFrame()");
439     GMX_RELEASE_ASSERT(!storedFrame.isFinished(),
440                        "finishFrame() called twice for the same frame");
441     GMX_RELEASE_ASSERT(storedFrame.frame->frameIndex() == index,
442                        "Inconsistent internal frame indexing");
443     storedFrame.status = Impl::StoredFrame::eFinished;
444     if (impl_->isMultipoint())
445     {
446         // TODO: Check that the last point set has been finished
447         impl_->data_->notifyFrameFinish(storedFrame.frame->header());
448         if (storedFrame.frame->frameIndex() >= impl_->storageLimit_)
449         {
450             impl_->rotateBuffer();
451         }
452     }
453     else
454     {
455         impl_->notifyNextFrames(storageIndex);
456     }
457 }
458
459
460 void
461 AnalysisDataStorage::finishFrame(const AnalysisDataStorageFrame &frame)
462 {
463     finishFrame(frame.frameIndex());
464 }
465
466 } // namespace gmx