2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2013,2014,2015,2016,2017,2018, 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.
37 * Implements classes in analysismodule.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_trajectoryanalysis
44 #include "analysismodule.h"
49 #include "gromacs/analysisdata/analysisdata.h"
50 #include "gromacs/selection/selection.h"
51 #include "gromacs/utility/exceptions.h"
52 #include "gromacs/utility/gmxassert.h"
57 /********************************************************************
58 * TrajectoryAnalysisModule::Impl
62 * Private implementation class for TrajectoryAnalysisModule.
64 * \ingroup module_trajectoryanalysis
66 class TrajectoryAnalysisModule::Impl
69 //! Container that associates a data set with its name.
70 typedef std::map
<std::string
, AbstractAnalysisData
*> DatasetContainer
;
71 //! Container that associates a AnalysisData object with its name.
72 typedef std::map
<std::string
, AnalysisData
*> AnalysisDatasetContainer
;
74 //! List of registered data set names.
75 std::vector
<std::string
> datasetNames_
;
77 * Keeps all registered data sets.
79 * This container also includes datasets from \a analysisDatasets_.
81 DatasetContainer datasets_
;
82 //! Keeps registered AnalysisData objects.
83 AnalysisDatasetContainer analysisDatasets_
;
86 /********************************************************************
87 * TrajectoryAnalysisModuleData::Impl
91 * Private implementation class for TrajectoryAnalysisModuleData.
93 * \ingroup module_trajectoryanalysis
95 class TrajectoryAnalysisModuleData::Impl
98 //! Container that associates a data handle to its AnalysisData object.
99 typedef std::map
<const AnalysisData
*, AnalysisDataHandle
>
102 //! \copydoc TrajectoryAnalysisModuleData::TrajectoryAnalysisModuleData()
103 Impl(TrajectoryAnalysisModule
*module
,
104 const AnalysisDataParallelOptions
&opt
,
105 const SelectionCollection
&selections
);
107 //! Checks whether the given AnalysisData has been initialized.
108 bool isInitialized(const AnalysisData
&data
) const;
110 //! Keeps a data handle for each AnalysisData object.
111 HandleContainer handles_
;
112 //! Stores thread-local selections.
113 const SelectionCollection
&selections_
;
116 TrajectoryAnalysisModuleData::Impl::Impl(
117 TrajectoryAnalysisModule
*module
,
118 const AnalysisDataParallelOptions
&opt
,
119 const SelectionCollection
&selections
)
120 : selections_(selections
)
122 TrajectoryAnalysisModule::Impl::AnalysisDatasetContainer::const_iterator i
;
123 for (i
= module
->impl_
->analysisDatasets_
.begin();
124 i
!= module
->impl_
->analysisDatasets_
.end(); ++i
)
126 AnalysisDataHandle handle
;
127 if (isInitialized(*i
->second
))
129 handle
= i
->second
->startData(opt
);
131 handles_
.insert(std::make_pair(i
->second
, handle
));
135 bool TrajectoryAnalysisModuleData::Impl::isInitialized(
136 const AnalysisData
&data
) const
138 for (int i
= 0; i
< data
.dataSetCount(); ++i
)
140 if (data
.columnCount(i
) > 0)
142 // If not all of the column counts are set, startData() in the
143 // constructor asserts, so that does not need to be checked here.
151 /********************************************************************
152 * TrajectoryAnalysisModuleData
155 TrajectoryAnalysisModuleData::TrajectoryAnalysisModuleData(
156 TrajectoryAnalysisModule
*module
,
157 const AnalysisDataParallelOptions
&opt
,
158 const SelectionCollection
&selections
)
159 : impl_(new Impl(module
, opt
, selections
))
164 TrajectoryAnalysisModuleData::~TrajectoryAnalysisModuleData()
169 void TrajectoryAnalysisModuleData::finishDataHandles()
171 // FIXME: Call finishData() for all handles even if one throws
172 Impl::HandleContainer::iterator i
;
173 for (i
= impl_
->handles_
.begin(); i
!= impl_
->handles_
.end(); ++i
)
175 if (i
->second
.isValid())
177 i
->second
.finishData();
180 impl_
->handles_
.clear();
185 TrajectoryAnalysisModuleData::dataHandle(const AnalysisData
&data
)
187 Impl::HandleContainer::const_iterator i
= impl_
->handles_
.find(&data
);
188 GMX_RELEASE_ASSERT(i
!= impl_
->handles_
.end(),
189 "Data handle requested on unknown dataset");
194 Selection
TrajectoryAnalysisModuleData::parallelSelection(const Selection
&selection
)
196 // TODO: Implement properly.
202 TrajectoryAnalysisModuleData::parallelSelections(const SelectionList
&selections
)
204 // TODO: Consider an implementation that does not allocate memory every time.
205 SelectionList newSelections
;
206 newSelections
.reserve(selections
.size());
207 SelectionList::const_iterator i
= selections
.begin();
208 for (; i
!= selections
.end(); ++i
)
210 newSelections
.push_back(parallelSelection(*i
));
212 return newSelections
;
216 /********************************************************************
217 * TrajectoryAnalysisModuleDataBasic
224 * Basic thread-local trajectory analysis data storage class.
226 * Most simple tools should only require data handles and selections to be
227 * thread-local, so this class implements just that.
229 * \ingroup module_trajectoryanalysis
231 class TrajectoryAnalysisModuleDataBasic
: public TrajectoryAnalysisModuleData
235 * Initializes thread-local storage for data handles and selections.
237 * \param[in] module Analysis module to use for data objects.
238 * \param[in] opt Data parallelization options.
239 * \param[in] selections Thread-local selection collection.
241 TrajectoryAnalysisModuleDataBasic(TrajectoryAnalysisModule
*module
,
242 const AnalysisDataParallelOptions
&opt
,
243 const SelectionCollection
&selections
);
245 void finish() override
;
248 TrajectoryAnalysisModuleDataBasic::TrajectoryAnalysisModuleDataBasic(
249 TrajectoryAnalysisModule
*module
,
250 const AnalysisDataParallelOptions
&opt
,
251 const SelectionCollection
&selections
)
252 : TrajectoryAnalysisModuleData(module
, opt
, selections
)
258 TrajectoryAnalysisModuleDataBasic::finish()
266 /********************************************************************
267 * TrajectoryAnalysisModule
270 TrajectoryAnalysisModule::TrajectoryAnalysisModule()
276 TrajectoryAnalysisModule::~TrajectoryAnalysisModule()
281 void TrajectoryAnalysisModule::optionsFinished(
282 TrajectoryAnalysisSettings
* /*settings*/)
287 void TrajectoryAnalysisModule::initAfterFirstFrame(
288 const TrajectoryAnalysisSettings
& /*settings*/,
289 const t_trxframe
& /*fr*/)
294 TrajectoryAnalysisModuleDataPointer
295 TrajectoryAnalysisModule::startFrames(const AnalysisDataParallelOptions
&opt
,
296 const SelectionCollection
&selections
)
298 return TrajectoryAnalysisModuleDataPointer(
299 new TrajectoryAnalysisModuleDataBasic(this, opt
, selections
));
303 void TrajectoryAnalysisModule::finishFrames(TrajectoryAnalysisModuleData
* /*pdata*/)
308 int TrajectoryAnalysisModule::datasetCount() const
310 return impl_
->datasetNames_
.size();
314 const std::vector
<std::string
> &TrajectoryAnalysisModule::datasetNames() const
316 return impl_
->datasetNames_
;
320 AbstractAnalysisData
&TrajectoryAnalysisModule::datasetFromIndex(int index
) const
322 if (index
< 0 || index
>= datasetCount())
324 GMX_THROW(APIError("Out of range data set index"));
326 Impl::DatasetContainer::const_iterator item
327 = impl_
->datasets_
.find(impl_
->datasetNames_
[index
]);
328 GMX_RELEASE_ASSERT(item
!= impl_
->datasets_
.end(),
329 "Inconsistent data set names");
330 return *item
->second
;
334 AbstractAnalysisData
&TrajectoryAnalysisModule::datasetFromName(const char *name
) const
336 Impl::DatasetContainer::const_iterator item
= impl_
->datasets_
.find(name
);
337 if (item
== impl_
->datasets_
.end())
339 GMX_THROW(APIError("Unknown data set name"));
341 return *item
->second
;
345 void TrajectoryAnalysisModule::registerBasicDataset(AbstractAnalysisData
*data
,
348 GMX_RELEASE_ASSERT(data
!= nullptr, "Attempting to register NULL data");
349 // TODO: Strong exception safety should be possible to implement.
350 GMX_RELEASE_ASSERT(impl_
->datasets_
.find(name
) == impl_
->datasets_
.end(),
351 "Duplicate data set name registered");
352 impl_
->datasets_
[name
] = data
;
353 impl_
->datasetNames_
.emplace_back(name
);
357 void TrajectoryAnalysisModule::registerAnalysisDataset(AnalysisData
*data
,
360 // TODO: Strong exception safety should be possible to implement.
361 registerBasicDataset(data
, name
);
362 impl_
->analysisDatasets_
[name
] = data
;
366 void TrajectoryAnalysisModule::finishFrameSerial(int frameIndex
)
368 Impl::AnalysisDatasetContainer::const_iterator data
;
369 for (data
= impl_
->analysisDatasets_
.begin();
370 data
!= impl_
->analysisDatasets_
.end();
373 data
->second
->finishFrameSerial(frameIndex
);