clang-tidy modernize
[gromacs.git] / src / gromacs / trajectoryanalysis / runnercommon.cpp
blob1b7e9a5ef75a6918e7ddcd5ef76c045978d78633
1 /*
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.
35 /*! \internal \file
36 * \brief
37 * Implements gmx::TrajectoryAnalysisRunnerCommon.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_trajectoryanalysis
42 #include "gmxpre.h"
44 #include "runnercommon.h"
46 #include <cstring>
48 #include <algorithm>
49 #include <string>
51 #include "gromacs/fileio/oenv.h"
52 #include "gromacs/fileio/timecontrol.h"
53 #include "gromacs/fileio/trxio.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/options/basicoptions.h"
56 #include "gromacs/options/filenameoption.h"
57 #include "gromacs/options/ioptionscontainer.h"
58 #include "gromacs/pbcutil/rmpbc.h"
59 #include "gromacs/selection/selection.h"
60 #include "gromacs/selection/selectioncollection.h"
61 #include "gromacs/selection/selectionoption.h"
62 #include "gromacs/selection/selectionoptionbehavior.h"
63 #include "gromacs/topology/topology.h"
64 #include "gromacs/trajectory/trajectoryframe.h"
65 #include "gromacs/trajectoryanalysis/analysissettings.h"
66 #include "gromacs/trajectoryanalysis/topologyinformation.h"
67 #include "gromacs/utility/cstringutil.h"
68 #include "gromacs/utility/exceptions.h"
69 #include "gromacs/utility/gmxassert.h"
70 #include "gromacs/utility/programcontext.h"
71 #include "gromacs/utility/smalloc.h"
72 #include "gromacs/utility/stringutil.h"
74 #include "analysissettings-impl.h"
76 namespace gmx
79 class TrajectoryAnalysisRunnerCommon::Impl : public ITopologyProvider
81 public:
82 explicit Impl(TrajectoryAnalysisSettings *settings);
83 ~Impl();
85 bool hasTrajectory() const { return !trjfile_.empty(); }
87 void initTopology(bool required);
88 void initFirstFrame();
89 void initFrameIndexGroup();
90 void finishTrajectory();
92 // From ITopologyProvider
93 virtual gmx_mtop_t *getTopology(bool required)
95 initTopology(required);
96 return topInfo_.mtop_.get();
98 virtual int getAtomCount()
100 if (!topInfo_.hasTopology())
102 if (trajectoryGroup_.isValid())
104 GMX_THROW(InconsistentInputError("-fgroup is only supported when -s is also specified"));
106 // Read the first frame if we don't know the maximum number of
107 // atoms otherwise.
108 initFirstFrame();
109 return fr->natoms;
111 return -1;
114 TrajectoryAnalysisSettings &settings_;
115 TopologyInformation topInfo_;
117 //! Name of the trajectory file (empty if not provided).
118 std::string trjfile_;
119 //! Name of the topology file (empty if no topology provided).
120 std::string topfile_;
121 Selection trajectoryGroup_;
122 double startTime_;
123 double endTime_;
124 double deltaTime_;
125 bool bStartTimeSet_;
126 bool bEndTimeSet_;
127 bool bDeltaTimeSet_;
129 bool bTrajOpen_;
130 //! The current frame, or \p NULL if no frame loaded yet.
131 t_trxframe *fr;
132 gmx_rmpbc_t gpbc_;
133 //! Used to store the status variable from read_first_frame().
134 t_trxstatus *status_;
135 gmx_output_env_t *oenv_;
139 TrajectoryAnalysisRunnerCommon::Impl::Impl(TrajectoryAnalysisSettings *settings)
140 : settings_(*settings),
141 startTime_(0.0), endTime_(0.0), deltaTime_(0.0),
142 bStartTimeSet_(false), bEndTimeSet_(false), bDeltaTimeSet_(false),
143 bTrajOpen_(false), fr(nullptr), gpbc_(nullptr), status_(nullptr), oenv_(nullptr)
148 TrajectoryAnalysisRunnerCommon::Impl::~Impl()
150 finishTrajectory();
151 if (fr != nullptr)
153 // There doesn't seem to be a function for freeing frame data
154 sfree(fr->x);
155 sfree(fr->v);
156 sfree(fr->f);
157 sfree(fr->index);
158 sfree(fr);
160 if (oenv_ != nullptr)
162 output_env_done(oenv_);
166 void
167 TrajectoryAnalysisRunnerCommon::Impl::initTopology(bool required)
169 // Return immediately if the topology has already been loaded.
170 if (topInfo_.hasTopology())
172 return;
175 if (required && topfile_.empty())
177 GMX_THROW(InconsistentInputError("No topology provided, but one is required for analysis"));
180 // Load the topology if requested.
181 if (!topfile_.empty())
183 topInfo_.fillFromInputFile(topfile_);
184 if (hasTrajectory()
185 && !settings_.hasFlag(TrajectoryAnalysisSettings::efUseTopX))
187 sfree(topInfo_.xtop_);
188 topInfo_.xtop_ = nullptr;
193 void
194 TrajectoryAnalysisRunnerCommon::Impl::initFirstFrame()
196 // Return if we have already initialized the trajectory.
197 if (fr != nullptr)
199 return;
201 time_unit_t time_unit
202 = static_cast<time_unit_t>(settings_.timeUnit() + 1); // NOLINT(misc-misplaced-widening-cast)
203 output_env_init(&oenv_, getProgramContext(), time_unit, FALSE, exvgNONE, 0);
205 int frflags = settings_.frflags();
206 frflags |= TRX_NEED_X;
208 snew(fr, 1);
210 if (hasTrajectory())
212 if (!read_first_frame(oenv_, &status_, trjfile_.c_str(), fr, frflags))
214 GMX_THROW(FileIOError("Could not read coordinates from trajectory"));
216 bTrajOpen_ = true;
218 if (topInfo_.hasTopology())
220 const int topologyAtomCount = topInfo_.topology()->atoms.nr;
221 if (fr->natoms > topologyAtomCount)
223 const std::string message
224 = formatString("Trajectory (%d atoms) does not match topology (%d atoms)",
225 fr->natoms, topologyAtomCount);
226 GMX_THROW(InconsistentInputError(message));
230 else
232 // Prepare a frame from topology information.
233 // TODO: Initialize more of the fields.
234 if (frflags & (TRX_NEED_V))
236 GMX_THROW(NotImplementedError("Velocity reading from a topology not implemented"));
238 if (frflags & (TRX_NEED_F))
240 GMX_THROW(InvalidInputError("Forces cannot be read from a topology"));
242 fr->natoms = topInfo_.topology()->atoms.nr;
243 fr->bX = TRUE;
244 snew(fr->x, fr->natoms);
245 memcpy(fr->x, topInfo_.xtop_,
246 sizeof(*fr->x) * fr->natoms);
247 fr->bBox = TRUE;
248 copy_mat(topInfo_.boxtop_, fr->box);
251 set_trxframe_ePBC(fr, topInfo_.ePBC());
252 if (topInfo_.hasTopology() && settings_.hasRmPBC())
254 gpbc_ = gmx_rmpbc_init(&topInfo_.topology()->idef, topInfo_.ePBC(),
255 fr->natoms);
259 void
260 TrajectoryAnalysisRunnerCommon::Impl::initFrameIndexGroup()
262 if (!trajectoryGroup_.isValid())
264 return;
266 GMX_RELEASE_ASSERT(bTrajOpen_,
267 "Trajectory index only makes sense with a real trajectory");
268 if (trajectoryGroup_.atomCount() != fr->natoms)
270 const std::string message = formatString(
271 "Selection specified with -fgroup has %d atoms, but "
272 "the trajectory (-f) has %d atoms.",
273 trajectoryGroup_.atomCount(), fr->natoms);
274 GMX_THROW(InconsistentInputError(message));
276 fr->bIndex = TRUE;
277 snew(fr->index, trajectoryGroup_.atomCount());
278 std::copy(trajectoryGroup_.atomIndices().begin(),
279 trajectoryGroup_.atomIndices().end(),
280 fr->index);
283 void
284 TrajectoryAnalysisRunnerCommon::Impl::finishTrajectory()
286 if (bTrajOpen_)
288 close_trx(status_);
289 bTrajOpen_ = false;
291 if (gpbc_ != nullptr)
293 gmx_rmpbc_done(gpbc_);
294 gpbc_ = nullptr;
298 /*********************************************************************
299 * TrajectoryAnalysisRunnerCommon
302 TrajectoryAnalysisRunnerCommon::TrajectoryAnalysisRunnerCommon(
303 TrajectoryAnalysisSettings *settings)
304 : impl_(new Impl(settings))
309 TrajectoryAnalysisRunnerCommon::~TrajectoryAnalysisRunnerCommon()
314 ITopologyProvider *
315 TrajectoryAnalysisRunnerCommon::topologyProvider()
317 return impl_.get();
321 void
322 TrajectoryAnalysisRunnerCommon::initOptions(IOptionsContainer *options,
323 TimeUnitBehavior *timeUnitBehavior)
325 TrajectoryAnalysisSettings &settings = impl_->settings_;
327 // Add common file name arguments.
328 options->addOption(FileNameOption("f")
329 .filetype(eftTrajectory).inputFile()
330 .store(&impl_->trjfile_)
331 .defaultBasename("traj")
332 .description("Input trajectory or single configuration"));
333 options->addOption(FileNameOption("s")
334 .filetype(eftTopology).inputFile()
335 .store(&impl_->topfile_)
336 .defaultBasename("topol")
337 .description("Input structure"));
339 // Add options for trajectory time control.
340 options->addOption(DoubleOption("b")
341 .store(&impl_->startTime_).storeIsSet(&impl_->bStartTimeSet_)
342 .timeValue()
343 .description("First frame (%t) to read from trajectory"));
344 options->addOption(DoubleOption("e")
345 .store(&impl_->endTime_).storeIsSet(&impl_->bEndTimeSet_)
346 .timeValue()
347 .description("Last frame (%t) to read from trajectory"));
348 options->addOption(DoubleOption("dt")
349 .store(&impl_->deltaTime_).storeIsSet(&impl_->bDeltaTimeSet_)
350 .timeValue()
351 .description("Only use frame if t MOD dt == first time (%t)"));
353 // Add time unit option.
354 timeUnitBehavior->setTimeUnitFromEnvironment();
355 timeUnitBehavior->addTimeUnitOption(options, "tu");
356 timeUnitBehavior->setTimeUnitStore(&impl_->settings_.impl_->timeUnit);
358 options->addOption(SelectionOption("fgroup")
359 .store(&impl_->trajectoryGroup_)
360 .onlySortedAtoms().onlyStatic()
361 .description("Atoms stored in the trajectory file "
362 "(if not set, assume first N atoms)"));
364 // Add plot options.
365 settings.impl_->plotSettings.initOptions(options);
367 // Add common options for trajectory processing.
368 if (!settings.hasFlag(TrajectoryAnalysisSettings::efNoUserRmPBC))
370 options->addOption(BooleanOption("rmpbc").store(&settings.impl_->bRmPBC)
371 .description("Make molecules whole for each frame"));
373 if (!settings.hasFlag(TrajectoryAnalysisSettings::efNoUserPBC))
375 options->addOption(BooleanOption("pbc").store(&settings.impl_->bPBC)
376 .description("Use periodic boundary conditions for distance calculation"));
381 void
382 TrajectoryAnalysisRunnerCommon::optionsFinished()
384 if (impl_->trjfile_.empty() && impl_->topfile_.empty())
386 GMX_THROW(InconsistentInputError("No trajectory or topology provided, nothing to do!"));
389 if (impl_->trajectoryGroup_.isValid() && impl_->trjfile_.empty())
391 GMX_THROW(InconsistentInputError("-fgroup only makes sense together with a trajectory (-f)"));
394 impl_->settings_.impl_->plotSettings.setTimeUnit(impl_->settings_.timeUnit());
396 if (impl_->bStartTimeSet_)
398 setTimeValue(TBEGIN, impl_->startTime_);
400 if (impl_->bEndTimeSet_)
402 setTimeValue(TEND, impl_->endTime_);
404 if (impl_->bDeltaTimeSet_)
406 setTimeValue(TDELTA, impl_->deltaTime_);
411 void
412 TrajectoryAnalysisRunnerCommon::initTopology()
414 const bool topologyRequired =
415 impl_->settings_.hasFlag(TrajectoryAnalysisSettings::efRequireTop);
416 impl_->initTopology(topologyRequired);
420 void
421 TrajectoryAnalysisRunnerCommon::initFirstFrame()
423 impl_->initFirstFrame();
427 void
428 TrajectoryAnalysisRunnerCommon::initFrameIndexGroup()
430 impl_->initFrameIndexGroup();
434 bool
435 TrajectoryAnalysisRunnerCommon::readNextFrame()
437 bool bContinue = false;
438 if (hasTrajectory())
440 bContinue = read_next_frame(impl_->oenv_, impl_->status_, impl_->fr);
442 if (!bContinue)
444 impl_->finishTrajectory();
446 return bContinue;
450 void
451 TrajectoryAnalysisRunnerCommon::initFrame()
453 if (impl_->gpbc_ != nullptr)
455 gmx_rmpbc_trxfr(impl_->gpbc_, impl_->fr);
460 bool
461 TrajectoryAnalysisRunnerCommon::hasTrajectory() const
463 return impl_->hasTrajectory();
467 const TopologyInformation &
468 TrajectoryAnalysisRunnerCommon::topologyInformation() const
470 return impl_->topInfo_;
474 t_trxframe &
475 TrajectoryAnalysisRunnerCommon::frame() const
477 GMX_RELEASE_ASSERT(impl_->fr != nullptr, "Frame not available when accessed");
478 return *impl_->fr;
481 } // namespace gmx