Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / analysisdata / modules / plot.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements classes in plot.h.
38  *
39  * \ingroup module_analysisdata
40  * \author Teemu Murtola <teemu.murtola@gmail.com>
41  */
42 #include "gmxpre.h"
43
44 #include "plot.h"
45
46 #include <cstdio>
47 #include <cstring>
48
49 #include <string>
50 #include <vector>
51
52 #include <boost/shared_ptr.hpp>
53
54 #include "gromacs/analysisdata/dataframe.h"
55 #include "gromacs/fileio/gmxfio.h"
56 #include "gromacs/fileio/xvgr.h"
57 #include "gromacs/legacyheaders/oenv.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/options/basicoptions.h"
60 #include "gromacs/options/options.h"
61 #include "gromacs/options/timeunitmanager.h"
62 #include "gromacs/selection/selectioncollection.h"
63 #include "gromacs/utility/exceptions.h"
64 #include "gromacs/utility/gmxassert.h"
65 #include "gromacs/utility/programcontext.h"
66 #include "gromacs/utility/stringutil.h"
67
68 namespace
69 {
70
71 //! Enum values for plot formats.
72 const char *const g_plotFormats[] = {
73     "none", "xmgrace", "xmgr"
74 };
75
76 } // namespace
77
78 namespace gmx
79 {
80
81 /********************************************************************
82  * AnalysisDataPlotSettings
83  */
84
85 AnalysisDataPlotSettings::AnalysisDataPlotSettings()
86     : selections_(NULL), timeUnit_(eTimeUnit_ps), plotFormat_(1)
87 {
88 }
89
90 void
91 AnalysisDataPlotSettings::setSelectionCollection(const SelectionCollection *selections)
92 {
93     selections_ = selections;
94 }
95
96
97 void
98 AnalysisDataPlotSettings::addOptions(Options *options)
99 {
100     options->addOption(StringOption("xvg").enumValue(g_plotFormats)
101                            .defaultValue("xmgrace")
102                            .storeEnumIndex(&plotFormat_)
103                            .description("Plot formatting"));
104 }
105
106
107 /********************************************************************
108  * AbstractPlotModule::Impl
109  */
110
111 class AbstractPlotModule::Impl
112 {
113     public:
114         explicit Impl(const AnalysisDataPlotSettings &settings);
115         ~Impl();
116
117         void closeFile();
118
119         AnalysisDataPlotSettings  settings_;
120         std::string               filename_;
121         FILE                     *fp_;
122
123         bool                      bPlain_;
124         bool                      bOmitX_;
125         bool                      bErrorsAsSeparateColumn_;
126         std::string               title_;
127         std::string               subtitle_;
128         std::string               xlabel_;
129         std::string               ylabel_;
130         std::vector<std::string>  legend_;
131         char                      xformat_[15];
132         char                      yformat_[15];
133         real                      xscale_;
134 };
135
136 AbstractPlotModule::Impl::Impl(const AnalysisDataPlotSettings &settings)
137     : settings_(settings), fp_(NULL), bPlain_(false), bOmitX_(false),
138       bErrorsAsSeparateColumn_(false), xscale_(1.0)
139 {
140     strcpy(xformat_, "%11.3f");
141     strcpy(yformat_, " %8.3f");
142 }
143
144 AbstractPlotModule::Impl::~Impl()
145 {
146     closeFile();
147 }
148
149
150 void
151 AbstractPlotModule::Impl::closeFile()
152 {
153     if (fp_ != NULL)
154     {
155         if (bPlain_)
156         {
157             gmx_fio_fclose(fp_);
158         }
159         else
160         {
161             xvgrclose(fp_);
162         }
163         fp_ = NULL;
164     }
165 }
166
167
168 /********************************************************************
169  * AbstractPlotModule
170  */
171 /*! \cond libapi */
172 AbstractPlotModule::AbstractPlotModule()
173     : impl_(new Impl(AnalysisDataPlotSettings()))
174 {
175 }
176
177 AbstractPlotModule::AbstractPlotModule(const AnalysisDataPlotSettings &settings)
178     : impl_(new Impl(settings))
179 {
180 }
181 //! \endcond
182
183 AbstractPlotModule::~AbstractPlotModule()
184 {
185 }
186
187
188 void
189 AbstractPlotModule::setSettings(const AnalysisDataPlotSettings &settings)
190 {
191     impl_->settings_ = settings;
192 }
193
194
195 void
196 AbstractPlotModule::setFileName(const std::string &filename)
197 {
198     impl_->filename_ = filename;
199 }
200
201
202 void
203 AbstractPlotModule::setPlainOutput(bool bPlain)
204 {
205     impl_->bPlain_ = bPlain;
206 }
207
208
209 void
210 AbstractPlotModule::setErrorsAsSeparateColumn(bool bSeparate)
211 {
212     impl_->bErrorsAsSeparateColumn_ = bSeparate;
213 }
214
215
216 void
217 AbstractPlotModule::setOmitX(bool bOmitX)
218 {
219     impl_->bOmitX_ = bOmitX;
220 }
221
222
223 void
224 AbstractPlotModule::setTitle(const char *title)
225 {
226     impl_->title_ = title;
227 }
228
229
230 void
231 AbstractPlotModule::setSubtitle(const char *subtitle)
232 {
233     impl_->subtitle_ = subtitle;
234 }
235
236
237 void
238 AbstractPlotModule::setXLabel(const char *label)
239 {
240     impl_->xlabel_ = label;
241 }
242
243
244 void
245 AbstractPlotModule::setXAxisIsTime()
246 {
247     TimeUnitManager manager(impl_->settings_.timeUnit());
248     impl_->xlabel_ = formatString("Time (%s)", manager.timeUnitAsString());
249     impl_->xscale_ = manager.inverseTimeScaleFactor();
250 }
251
252
253 void
254 AbstractPlotModule::setYLabel(const char *label)
255 {
256     impl_->ylabel_ = label;
257 }
258
259
260 void
261 AbstractPlotModule::setLegend(int nsets, const char * const *setname)
262 {
263     impl_->legend_.reserve(impl_->legend_.size() + nsets);
264     for (int i = 0; i < nsets; ++i)
265     {
266         appendLegend(setname[i]);
267     }
268 }
269
270
271 void
272 AbstractPlotModule::appendLegend(const char *setname)
273 {
274     impl_->legend_.push_back(setname);
275 }
276
277
278 void
279 AbstractPlotModule::appendLegend(const std::string &setname)
280 {
281     impl_->legend_.push_back(setname);
282 }
283
284
285 void
286 AbstractPlotModule::setXFormat(int width, int precision, char format)
287 {
288     GMX_RELEASE_ASSERT(width >= 0 && precision >= 0
289                        && width <= 99 && precision <= 99,
290                        "Invalid width or precision");
291     GMX_RELEASE_ASSERT(strchr("eEfFgG", format) != NULL,
292                        "Invalid format specifier");
293     sprintf(impl_->xformat_, "%%%d.%d%c", width, precision, format);
294 }
295
296
297 void
298 AbstractPlotModule::setYFormat(int width, int precision, char format)
299 {
300     GMX_RELEASE_ASSERT(width >= 0 && precision >= 0
301                        && width <= 99 && precision <= 99,
302                        "Invalid width or precision");
303     GMX_RELEASE_ASSERT(strchr("eEfFgG", format) != NULL,
304                        "Invalid format specifier");
305     sprintf(impl_->yformat_, " %%%d.%d%c", width, precision, format);
306 }
307
308
309 int
310 AbstractPlotModule::flags() const
311 {
312     return efAllowMissing | efAllowMulticolumn | efAllowMultipoint
313            | efAllowMultipleDataSets;
314 }
315
316
317 void
318 AbstractPlotModule::dataStarted(AbstractAnalysisData * /* data */)
319 {
320     if (!impl_->filename_.empty())
321     {
322         if (impl_->bPlain_)
323         {
324             impl_->fp_ = gmx_fio_fopen(impl_->filename_.c_str(), "w");
325         }
326         else
327         {
328             time_unit_t  time_unit
329                 = static_cast<time_unit_t>(impl_->settings_.timeUnit() + 1);
330             xvg_format_t xvg_format
331                 = (impl_->settings_.plotFormat() > 0
332                    ? static_cast<xvg_format_t>(impl_->settings_.plotFormat())
333                    : exvgNONE);
334             output_env_t                  oenv;
335             output_env_init(&oenv, getProgramContext(), time_unit, FALSE, xvg_format, 0);
336             boost::shared_ptr<output_env> oenvGuard(oenv, &output_env_done);
337             impl_->fp_ = xvgropen(impl_->filename_.c_str(), impl_->title_.c_str(),
338                                   impl_->xlabel_.c_str(), impl_->ylabel_.c_str(),
339                                   oenv);
340             const SelectionCollection *selections
341                 = impl_->settings_.selectionCollection();
342             if (selections != NULL)
343             {
344                 selections->printXvgrInfo(impl_->fp_, oenv);
345             }
346             if (!impl_->subtitle_.empty())
347             {
348                 xvgr_subtitle(impl_->fp_, impl_->subtitle_.c_str(), oenv);
349             }
350             if (output_env_get_print_xvgr_codes(oenv)
351                 && !impl_->legend_.empty())
352             {
353                 std::vector<const char *> legend;
354                 legend.reserve(impl_->legend_.size());
355                 for (size_t i = 0; i < impl_->legend_.size(); ++i)
356                 {
357                     legend.push_back(impl_->legend_[i].c_str());
358                 }
359                 xvgr_legend(impl_->fp_, legend.size(), &legend[0], oenv);
360             }
361         }
362     }
363 }
364
365
366 void
367 AbstractPlotModule::frameStarted(const AnalysisDataFrameHeader &frame)
368 {
369     if (!isFileOpen())
370     {
371         return;
372     }
373     if (!impl_->bOmitX_)
374     {
375         std::fprintf(impl_->fp_, impl_->xformat_, frame.x() * impl_->xscale_);
376     }
377 }
378
379
380 void
381 AbstractPlotModule::frameFinished(const AnalysisDataFrameHeader & /*header*/)
382 {
383     if (!isFileOpen())
384     {
385         return;
386     }
387     std::fprintf(impl_->fp_, "\n");
388 }
389
390
391 void
392 AbstractPlotModule::dataFinished()
393 {
394     impl_->closeFile();
395 }
396
397 /*! \cond libapi */
398 bool
399 AbstractPlotModule::isFileOpen() const
400 {
401     return impl_->fp_ != NULL;
402 }
403
404
405 void
406 AbstractPlotModule::writeValue(const AnalysisDataValue &value) const
407 {
408     GMX_ASSERT(isFileOpen(), "File not opened, but write attempted");
409     const real y = value.isSet() ? value.value() : 0.0;
410     std::fprintf(impl_->fp_, impl_->yformat_, y);
411     if (impl_->bErrorsAsSeparateColumn_)
412     {
413         const real dy = value.isSet() ? value.error() : 0.0;
414         std::fprintf(impl_->fp_, impl_->yformat_, dy);
415     }
416 }
417 //! \endcond
418
419 /********************************************************************
420  * DataPlotModule
421  */
422
423 AnalysisDataPlotModule::AnalysisDataPlotModule()
424 {
425 }
426
427 AnalysisDataPlotModule::AnalysisDataPlotModule(
428         const AnalysisDataPlotSettings &settings)
429     : AbstractPlotModule(settings)
430 {
431 }
432
433
434 void
435 AnalysisDataPlotModule::pointsAdded(const AnalysisDataPointSetRef &points)
436 {
437     if (!isFileOpen())
438     {
439         return;
440     }
441     for (int i = 0; i < points.columnCount(); ++i)
442     {
443         writeValue(points.values()[i]);
444     }
445 }
446
447
448 /********************************************************************
449  * DataVectorPlotModule
450  */
451
452 AnalysisDataVectorPlotModule::AnalysisDataVectorPlotModule()
453 {
454     for (int i = 0; i < DIM; ++i)
455     {
456         bWrite_[i] = true;
457     }
458     bWrite_[DIM] = false;
459 }
460
461
462 AnalysisDataVectorPlotModule::AnalysisDataVectorPlotModule(
463         const AnalysisDataPlotSettings &settings)
464     : AbstractPlotModule(settings)
465 {
466     for (int i = 0; i < DIM; ++i)
467     {
468         bWrite_[i] = true;
469     }
470     bWrite_[DIM] = false;
471 }
472
473
474 void
475 AnalysisDataVectorPlotModule::setWriteX(bool bWrite)
476 {
477     bWrite_[XX] = bWrite;
478 }
479
480
481 void
482 AnalysisDataVectorPlotModule::setWriteY(bool bWrite)
483 {
484     bWrite_[YY] = bWrite;
485 }
486
487
488 void
489 AnalysisDataVectorPlotModule::setWriteZ(bool bWrite)
490 {
491     bWrite_[ZZ] = bWrite;
492 }
493
494
495 void
496 AnalysisDataVectorPlotModule::setWriteNorm(bool bWrite)
497 {
498     bWrite_[DIM] = bWrite;
499 }
500
501
502 void
503 AnalysisDataVectorPlotModule::setWriteMask(bool bWrite[DIM + 1])
504 {
505     for (int i = 0; i < DIM + 1; ++i)
506     {
507         bWrite_[i] = bWrite[i];
508     }
509 }
510
511
512 void
513 AnalysisDataVectorPlotModule::pointsAdded(const AnalysisDataPointSetRef &points)
514 {
515     if (points.firstColumn() % DIM != 0 || points.columnCount() % DIM != 0)
516     {
517         GMX_THROW(APIError("Partial data points"));
518     }
519     if (!isFileOpen())
520     {
521         return;
522     }
523     for (int i = 0; i < points.columnCount(); i += 3)
524     {
525         for (int d = 0; d < DIM; ++d)
526         {
527             if (bWrite_[i])
528             {
529                 writeValue(points.values()[i + d]);
530             }
531         }
532         if (bWrite_[DIM])
533         {
534             const rvec        y = { points.y(i), points.y(i + 1), points.y(i + 2) };
535             AnalysisDataValue value(norm(y));
536             writeValue(value);
537         }
538     }
539 }
540
541 } // namespace gmx