Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / analysisdata / analysisdata.cpp
blob9ad7136a43893b73beb42ebe0fbc011fc485ab2c
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief
37 * Implements classes in analysisdata.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_analysisdata
42 #include "gmxpre.h"
44 #include "analysisdata.h"
46 #include <memory>
48 #include "gromacs/analysisdata/dataframe.h"
49 #include "gromacs/analysisdata/datastorage.h"
50 #include "gromacs/analysisdata/paralleloptions.h"
51 #include "gromacs/utility/exceptions.h"
52 #include "gromacs/utility/gmxassert.h"
54 namespace gmx
57 /********************************************************************
58 * AnalysisDataHandleImpl
61 namespace internal
64 /*! \internal \brief
65 * Private implementation class for AnalysisDataHandle.
67 * \ingroup module_analysisdata
69 class AnalysisDataHandleImpl
71 public:
72 //! Creates a handle associated with the given data object.
73 explicit AnalysisDataHandleImpl(AnalysisData *data)
74 : data_(*data), currentFrame_(nullptr)
78 //! The data object that this handle belongs to.
79 AnalysisData &data_;
80 //! Current storage frame object, or NULL if no current frame.
81 AnalysisDataStorageFrame *currentFrame_;
84 } // namespace internal
86 /********************************************************************
87 * AnalysisData::Impl
90 /*! \internal \brief
91 * Private implementation class for AnalysisData.
93 * \ingroup module_analysisdata
95 class AnalysisData::Impl
97 public:
98 //! Smart pointer type to manage a data handle implementation.
99 typedef std::unique_ptr<internal::AnalysisDataHandleImpl>
100 HandlePointer;
101 //! Shorthand for a list of data handles.
102 typedef std::vector<HandlePointer> HandleList;
104 //! Storage implementation.
105 AnalysisDataStorage storage_;
106 /*! \brief
107 * List of handles for this data object.
109 * Note that AnalysisDataHandle objects also contain (raw) pointers
110 * to these objects.
112 HandleList handles_;
115 /********************************************************************
116 * AnalysisData
119 AnalysisData::AnalysisData()
120 : impl_(new Impl)
125 AnalysisData::~AnalysisData()
130 void
131 AnalysisData::setDataSetCount(int dataSetCount)
133 GMX_RELEASE_ASSERT(impl_->handles_.empty(),
134 "Cannot change data dimensionality after creating handles");
135 AbstractAnalysisData::setDataSetCount(dataSetCount);
139 void
140 AnalysisData::setColumnCount(int dataSet, int columnCount)
142 GMX_RELEASE_ASSERT(impl_->handles_.empty(),
143 "Cannot change data dimensionality after creating handles");
144 AbstractAnalysisData::setColumnCount(dataSet, columnCount);
148 void
149 AnalysisData::setMultipoint(bool bMultipoint)
151 GMX_RELEASE_ASSERT(impl_->handles_.empty(),
152 "Cannot change data type after creating handles");
153 AbstractAnalysisData::setMultipoint(bMultipoint);
158 AnalysisData::frameCount() const
160 return impl_->storage_.frameCount();
164 AnalysisDataHandle
165 AnalysisData::startData(const AnalysisDataParallelOptions &opt)
167 GMX_RELEASE_ASSERT(impl_->handles_.size() < static_cast<unsigned>(opt.parallelizationFactor()),
168 "Too many calls to startData() compared to provided options");
169 if (impl_->handles_.empty())
171 impl_->storage_.startParallelDataStorage(this, &moduleManager(), opt);
174 Impl::HandlePointer handle(new internal::AnalysisDataHandleImpl(this));
175 impl_->handles_.push_back(std::move(handle));
176 return AnalysisDataHandle(impl_->handles_.back().get());
180 void
181 AnalysisData::finishFrameSerial(int frameIndex)
183 impl_->storage_.finishFrameSerial(frameIndex);
187 void
188 AnalysisData::finishData(AnalysisDataHandle handle)
190 Impl::HandleList::iterator i;
192 for (i = impl_->handles_.begin(); i != impl_->handles_.end(); ++i)
194 if (i->get() == handle.impl_)
196 break;
199 GMX_RELEASE_ASSERT(i != impl_->handles_.end(),
200 "finishData() called for an unknown handle");
202 impl_->handles_.erase(i);
204 if (impl_->handles_.empty())
206 impl_->storage_.finishDataStorage();
211 AnalysisDataFrameRef
212 AnalysisData::tryGetDataFrameInternal(int index) const
214 return impl_->storage_.tryGetDataFrame(index);
218 bool
219 AnalysisData::requestStorageInternal(int nframes)
221 return impl_->storage_.requestStorage(nframes);
225 /********************************************************************
226 * AnalysisDataHandle
229 AnalysisDataHandle::AnalysisDataHandle()
230 : impl_(nullptr)
235 AnalysisDataHandle::AnalysisDataHandle(internal::AnalysisDataHandleImpl *impl)
236 : impl_(impl)
241 void
242 AnalysisDataHandle::startFrame(int index, real x, real dx)
244 GMX_RELEASE_ASSERT(impl_ != nullptr, "Invalid data handle used");
245 GMX_RELEASE_ASSERT(impl_->currentFrame_ == nullptr,
246 "startFrame() called twice without calling finishFrame()");
247 impl_->currentFrame_ =
248 &impl_->data_.impl_->storage_.startFrame(index, x, dx);
252 void
253 AnalysisDataHandle::selectDataSet(int index)
255 GMX_RELEASE_ASSERT(impl_ != nullptr, "Invalid data handle used");
256 GMX_RELEASE_ASSERT(impl_->currentFrame_ != nullptr,
257 "selectDataSet() called without calling startFrame()");
258 impl_->currentFrame_->selectDataSet(index);
262 void
263 AnalysisDataHandle::setPoint(int column, real value, bool bPresent)
265 GMX_RELEASE_ASSERT(impl_ != nullptr, "Invalid data handle used");
266 GMX_RELEASE_ASSERT(impl_->currentFrame_ != nullptr,
267 "setPoint() called without calling startFrame()");
268 impl_->currentFrame_->setValue(column, value, bPresent);
272 void
273 AnalysisDataHandle::setPoint(int column, real value, real error, bool bPresent)
275 GMX_RELEASE_ASSERT(impl_ != nullptr, "Invalid data handle used");
276 GMX_RELEASE_ASSERT(impl_->currentFrame_ != nullptr,
277 "setPoint() called without calling startFrame()");
278 impl_->currentFrame_->setValue(column, value, error, bPresent);
282 void
283 AnalysisDataHandle::setPoints(int firstColumn, int count, const real *values,
284 bool bPresent)
286 GMX_RELEASE_ASSERT(impl_ != nullptr, "Invalid data handle used");
287 GMX_RELEASE_ASSERT(impl_->currentFrame_ != nullptr,
288 "setPoints() called without calling startFrame()");
289 for (int i = 0; i < count; ++i)
291 impl_->currentFrame_->setValue(firstColumn + i, values[i], bPresent);
296 void
297 AnalysisDataHandle::finishPointSet()
299 GMX_RELEASE_ASSERT(impl_ != nullptr, "Invalid data handle used");
300 GMX_RELEASE_ASSERT(impl_->data_.isMultipoint(),
301 "finishPointSet() called for non-multipoint data");
302 GMX_RELEASE_ASSERT(impl_->currentFrame_ != nullptr,
303 "finishPointSet() called without calling startFrame()");
304 impl_->currentFrame_->finishPointSet();
308 void
309 AnalysisDataHandle::finishFrame()
311 GMX_RELEASE_ASSERT(impl_ != nullptr, "Invalid data handle used");
312 GMX_RELEASE_ASSERT(impl_->currentFrame_ != nullptr,
313 "finishFrame() called without calling startFrame()");
314 AnalysisDataStorageFrame *frame = impl_->currentFrame_;
315 impl_->currentFrame_ = nullptr;
316 frame->finishFrame();
320 void
321 AnalysisDataHandle::finishData()
323 GMX_RELEASE_ASSERT(impl_ != nullptr, "Invalid data handle used");
324 // Deletes the implementation pointer.
325 impl_->data_.finishData(*this);
326 impl_ = nullptr;
329 } // namespace gmx