Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / analysisdata / modules / plot.cpp
blob25eccc9bd4d06b1ad11ee8bdfc6c15575510f698
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 classes in plot.h.
39 * \ingroup module_analysisdata
40 * \author Teemu Murtola <teemu.murtola@gmail.com>
42 #include "gmxpre.h"
44 #include "plot.h"
46 #include <cstdio>
47 #include <cstring>
49 #include <string>
50 #include <vector>
52 #include "gromacs/analysisdata/dataframe.h"
53 #include "gromacs/fileio/gmxfio.h"
54 #include "gromacs/fileio/oenv.h"
55 #include "gromacs/fileio/xvgr.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/options/basicoptions.h"
58 #include "gromacs/options/ioptionscontainer.h"
59 #include "gromacs/options/timeunitmanager.h"
60 #include "gromacs/selection/selectioncollection.h"
61 #include "gromacs/utility/exceptions.h"
62 #include "gromacs/utility/gmxassert.h"
63 #include "gromacs/utility/programcontext.h"
64 #include "gromacs/utility/stringutil.h"
65 #include "gromacs/utility/unique_cptr.h"
67 namespace
70 //! Enum values for plot formats.
71 const char *const g_plotFormats[] = {
72 "none", "xmgrace", "xmgr"
75 } // namespace
77 namespace gmx
80 /********************************************************************
81 * AnalysisDataPlotSettings
84 AnalysisDataPlotSettings::AnalysisDataPlotSettings()
85 : selections_(nullptr), timeUnit_(TimeUnit_Default), plotFormat_(1)
89 void
90 AnalysisDataPlotSettings::setSelectionCollection(const SelectionCollection *selections)
92 selections_ = selections;
96 void
97 AnalysisDataPlotSettings::initOptions(IOptionsContainer *options)
99 options->addOption(EnumIntOption("xvg").enumValue(g_plotFormats)
100 .store(&plotFormat_)
101 .description("Plot formatting"));
105 /********************************************************************
106 * AbstractPlotModule::Impl
109 class AbstractPlotModule::Impl
111 public:
112 explicit Impl(const AnalysisDataPlotSettings &settings);
113 ~Impl();
115 void closeFile();
117 AnalysisDataPlotSettings settings_;
118 std::string filename_;
119 FILE *fp_;
121 bool bPlain_;
122 bool bOmitX_;
123 bool bErrorsAsSeparateColumn_;
124 std::string title_;
125 std::string subtitle_;
126 std::string xlabel_;
127 std::string ylabel_;
128 std::vector<std::string> legend_;
129 std::string xformat_;
130 std::string yformat_;
131 real xscale_;
134 AbstractPlotModule::Impl::Impl(const AnalysisDataPlotSettings &settings)
135 : settings_(settings), fp_(nullptr), bPlain_(false), bOmitX_(false),
136 bErrorsAsSeparateColumn_(false),
137 xformat_("%11.3f"), yformat_(" %8.3f"), xscale_(1.0)
141 AbstractPlotModule::Impl::~Impl()
143 closeFile();
147 void
148 AbstractPlotModule::Impl::closeFile()
150 if (fp_ != nullptr)
152 if (bPlain_)
154 gmx_fio_fclose(fp_);
156 else
158 xvgrclose(fp_);
160 fp_ = nullptr;
165 /********************************************************************
166 * AbstractPlotModule
168 /*! \cond libapi */
169 AbstractPlotModule::AbstractPlotModule()
170 : impl_(new Impl(AnalysisDataPlotSettings()))
174 AbstractPlotModule::AbstractPlotModule(const AnalysisDataPlotSettings &settings)
175 : impl_(new Impl(settings))
178 //! \endcond
180 AbstractPlotModule::~AbstractPlotModule()
185 void
186 AbstractPlotModule::setSettings(const AnalysisDataPlotSettings &settings)
188 impl_->settings_ = settings;
192 void
193 AbstractPlotModule::setFileName(const std::string &filename)
195 impl_->filename_ = filename;
199 void
200 AbstractPlotModule::setPlainOutput(bool bPlain)
202 impl_->bPlain_ = bPlain;
206 void
207 AbstractPlotModule::setErrorsAsSeparateColumn(bool bSeparate)
209 impl_->bErrorsAsSeparateColumn_ = bSeparate;
213 void
214 AbstractPlotModule::setOmitX(bool bOmitX)
216 impl_->bOmitX_ = bOmitX;
220 void
221 AbstractPlotModule::setTitle(const char *title)
223 impl_->title_ = title;
226 void
227 AbstractPlotModule::setTitle(const std::string &title)
229 impl_->title_ = title;
233 void
234 AbstractPlotModule::setSubtitle(const char *subtitle)
236 impl_->subtitle_ = subtitle;
240 void
241 AbstractPlotModule::setSubtitle(const std::string &subtitle)
243 impl_->subtitle_ = subtitle;
247 void
248 AbstractPlotModule::setXLabel(const char *label)
250 impl_->xlabel_ = label;
254 void
255 AbstractPlotModule::setXAxisIsTime()
257 TimeUnitManager manager(impl_->settings_.timeUnit());
258 impl_->xlabel_ = formatString("Time (%s)", manager.timeUnitAsString());
259 impl_->xscale_ = manager.inverseTimeScaleFactor();
263 void
264 AbstractPlotModule::setYLabel(const char *label)
266 impl_->ylabel_ = label;
270 void
271 AbstractPlotModule::setLegend(int nsets, const char * const *setname)
273 impl_->legend_.reserve(impl_->legend_.size() + nsets);
274 for (int i = 0; i < nsets; ++i)
276 appendLegend(setname[i]);
281 void
282 AbstractPlotModule::appendLegend(const char *setname)
284 impl_->legend_.emplace_back(setname);
288 void
289 AbstractPlotModule::appendLegend(const std::string &setname)
291 impl_->legend_.push_back(setname);
295 void
296 AbstractPlotModule::setXFormat(int width, int precision, char format)
298 GMX_RELEASE_ASSERT(width >= 0 && precision >= 0
299 && width <= 99 && precision <= 99,
300 "Invalid width or precision");
301 GMX_RELEASE_ASSERT(strchr("eEfFgG", format) != nullptr,
302 "Invalid format specifier");
303 impl_->xformat_ = formatString("%%%d.%d%c", width, precision, format);
307 void
308 AbstractPlotModule::setYFormat(int width, int precision, char format)
310 GMX_RELEASE_ASSERT(width >= 0 && precision >= 0
311 && width <= 99 && precision <= 99,
312 "Invalid width or precision");
313 GMX_RELEASE_ASSERT(strchr("eEfFgG", format) != nullptr,
314 "Invalid format specifier");
315 impl_->yformat_ = formatString(" %%%d.%d%c", width, precision, format);
320 AbstractPlotModule::flags() const
322 return efAllowMissing | efAllowMulticolumn | efAllowMultipoint
323 | efAllowMultipleDataSets;
327 void
328 AbstractPlotModule::dataStarted(AbstractAnalysisData * /* data */)
330 if (!impl_->filename_.empty())
332 if (impl_->bPlain_)
334 impl_->fp_ = gmx_fio_fopen(impl_->filename_.c_str(), "w");
336 else
338 time_unit_t time_unit
339 = static_cast<time_unit_t>(impl_->settings_.timeUnit() + 1); // NOLINT(bugprone-misplaced-widening-cast)
340 xvg_format_t xvg_format
341 = (impl_->settings_.plotFormat() > 0
342 ? static_cast<xvg_format_t>(impl_->settings_.plotFormat())
343 : exvgNONE);
344 gmx_output_env_t *oenv;
345 output_env_init(&oenv, getProgramContext(), time_unit, FALSE, xvg_format, 0);
346 const unique_cptr<gmx_output_env_t, output_env_done> oenvGuard(oenv);
347 impl_->fp_ = xvgropen(impl_->filename_.c_str(), impl_->title_.c_str(),
348 impl_->xlabel_, impl_->ylabel_,
349 oenv);
350 const SelectionCollection *selections
351 = impl_->settings_.selectionCollection();
352 if (selections != nullptr && output_env_get_xvg_format(oenv) != exvgNONE)
354 selections->printXvgrInfo(impl_->fp_);
356 if (!impl_->subtitle_.empty())
358 xvgr_subtitle(impl_->fp_, impl_->subtitle_.c_str(), oenv);
360 if (output_env_get_print_xvgr_codes(oenv)
361 && !impl_->legend_.empty())
363 std::vector<const char *> legend;
364 legend.reserve(impl_->legend_.size());
365 for (size_t i = 0; i < impl_->legend_.size(); ++i)
367 legend.push_back(impl_->legend_[i].c_str());
369 xvgr_legend(impl_->fp_, legend.size(), legend.data(), oenv);
376 void
377 AbstractPlotModule::frameStarted(const AnalysisDataFrameHeader &header)
379 if (!isFileOpen())
381 return;
383 if (!impl_->bOmitX_)
385 std::fprintf(impl_->fp_, impl_->xformat_.c_str(), header.x() * impl_->xscale_);
390 void
391 AbstractPlotModule::frameFinished(const AnalysisDataFrameHeader & /*header*/)
393 if (!isFileOpen())
395 return;
397 std::fprintf(impl_->fp_, "\n");
401 void
402 AbstractPlotModule::dataFinished()
404 impl_->closeFile();
407 /*! \cond libapi */
408 bool
409 AbstractPlotModule::isFileOpen() const
411 return impl_->fp_ != nullptr;
415 void
416 AbstractPlotModule::writeValue(const AnalysisDataValue &value) const
418 GMX_ASSERT(isFileOpen(), "File not opened, but write attempted");
419 const real y = value.isSet() ? value.value() : 0.0;
420 std::fprintf(impl_->fp_, impl_->yformat_.c_str(), y);
421 if (impl_->bErrorsAsSeparateColumn_)
423 const real dy = value.isSet() ? value.error() : 0.0;
424 std::fprintf(impl_->fp_, impl_->yformat_.c_str(), dy);
427 //! \endcond
429 /********************************************************************
430 * DataPlotModule
433 AnalysisDataPlotModule::AnalysisDataPlotModule()
437 AnalysisDataPlotModule::AnalysisDataPlotModule(
438 const AnalysisDataPlotSettings &settings)
439 : AbstractPlotModule(settings)
444 void
445 AnalysisDataPlotModule::pointsAdded(const AnalysisDataPointSetRef &points)
447 if (!isFileOpen())
449 return;
451 for (int i = 0; i < points.columnCount(); ++i)
453 writeValue(points.values()[i]);
458 /********************************************************************
459 * DataVectorPlotModule
462 AnalysisDataVectorPlotModule::AnalysisDataVectorPlotModule()
463 : bWrite_ {true, true, true, false}
468 AnalysisDataVectorPlotModule::AnalysisDataVectorPlotModule(
469 const AnalysisDataPlotSettings &settings)
470 : AbstractPlotModule(settings),
471 bWrite_ {true, true, true, false}
476 void
477 AnalysisDataVectorPlotModule::setWriteX(bool bWrite)
479 bWrite_[XX] = bWrite;
483 void
484 AnalysisDataVectorPlotModule::setWriteY(bool bWrite)
486 bWrite_[YY] = bWrite;
490 void
491 AnalysisDataVectorPlotModule::setWriteZ(bool bWrite)
493 bWrite_[ZZ] = bWrite;
497 void
498 AnalysisDataVectorPlotModule::setWriteNorm(bool bWrite)
500 bWrite_[DIM] = bWrite;
504 void
505 AnalysisDataVectorPlotModule::setWriteMask(const bool bWrite[DIM + 1])
507 for (int i = 0; i < DIM + 1; ++i)
509 bWrite_[i] = bWrite[i];
514 void
515 AnalysisDataVectorPlotModule::pointsAdded(const AnalysisDataPointSetRef &points)
517 if (points.firstColumn() % DIM != 0 || points.columnCount() % DIM != 0)
519 GMX_THROW(APIError("Partial data points"));
521 if (!isFileOpen())
523 return;
525 for (int i = 0; i < points.columnCount(); i += 3)
527 for (int d = 0; d < DIM; ++d)
529 if (bWrite_[d])
531 writeValue(points.values()[i + d]);
534 if (bWrite_[DIM])
536 const rvec y = { points.y(i), points.y(i + 1), points.y(i + 2) };
537 AnalysisDataValue value(norm(y));
538 writeValue(value);
543 } // namespace gmx