Version 5.1.5
[gromacs.git] / share / template / template.cpp
blob6ef9e466b5764566615f7be8b006efbea087ec71
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2011,2012,2013,2014,2015, 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 #include <string>
36 #include <vector>
38 #include <gromacs/trajectoryanalysis.h>
40 using namespace gmx;
42 /*! \brief
43 * Template class to serve as a basis for user analysis tools.
45 class AnalysisTemplate : public TrajectoryAnalysisModule
47 public:
48 AnalysisTemplate();
50 virtual void initOptions(Options *options,
51 TrajectoryAnalysisSettings *settings);
52 virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
53 const TopologyInformation &top);
55 virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
56 TrajectoryAnalysisModuleData *pdata);
58 virtual void finishAnalysis(int nframes);
59 virtual void writeOutput();
61 private:
62 class ModuleData;
64 std::string fnDist_;
65 double cutoff_;
66 Selection refsel_;
67 SelectionList sel_;
69 AnalysisNeighborhood nb_;
71 AnalysisData data_;
72 AnalysisDataAverageModulePointer avem_;
76 AnalysisTemplate::AnalysisTemplate()
77 : TrajectoryAnalysisModule("template", "Template analysis tool"),
78 cutoff_(0.0)
80 registerAnalysisDataset(&data_, "avedist");
84 void
85 AnalysisTemplate::initOptions(Options *options,
86 TrajectoryAnalysisSettings *settings)
88 static const char *const desc[] = {
89 "This is a template for writing your own analysis tools for",
90 "GROMACS. The advantage of using GROMACS for this is that you",
91 "have access to all information in the topology, and your",
92 "program will be able to handle all types of coordinates and",
93 "trajectory files supported by GROMACS. In addition,",
94 "you get a lot of functionality for free from the trajectory",
95 "analysis library, including support for flexible dynamic",
96 "selections. Go ahead an try it![PAR]",
97 "To get started with implementing your own analysis program,",
98 "follow the instructions in the README file provided.",
99 "This template implements a simple analysis programs that calculates",
100 "average distances from a reference group to one or more",
101 "analysis groups."
104 options->setDescription(desc);
106 options->addOption(FileNameOption("o")
107 .filetype(eftPlot).outputFile()
108 .store(&fnDist_).defaultBasename("avedist")
109 .description("Average distances from reference group"));
111 options->addOption(SelectionOption("reference")
112 .store(&refsel_).required()
113 .description("Reference group to calculate distances from"));
114 options->addOption(SelectionOption("select")
115 .storeVector(&sel_).required().multiValue()
116 .description("Groups to calculate distances to"));
118 options->addOption(DoubleOption("cutoff").store(&cutoff_)
119 .description("Cutoff for distance calculation (0 = no cutoff)"));
121 settings->setFlag(TrajectoryAnalysisSettings::efRequireTop);
125 void
126 AnalysisTemplate::initAnalysis(const TrajectoryAnalysisSettings &settings,
127 const TopologyInformation & /*top*/)
129 nb_.setCutoff(cutoff_);
131 data_.setColumnCount(0, sel_.size());
133 avem_.reset(new AnalysisDataAverageModule());
134 data_.addModule(avem_);
136 if (!fnDist_.empty())
138 AnalysisDataPlotModulePointer plotm(
139 new AnalysisDataPlotModule(settings.plotSettings()));
140 plotm->setFileName(fnDist_);
141 plotm->setTitle("Average distance");
142 plotm->setXAxisIsTime();
143 plotm->setYLabel("Distance (nm)");
144 data_.addModule(plotm);
149 void
150 AnalysisTemplate::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
151 TrajectoryAnalysisModuleData *pdata)
153 AnalysisDataHandle dh = pdata->dataHandle(data_);
154 const Selection &refsel = pdata->parallelSelection(refsel_);
156 AnalysisNeighborhoodSearch nbsearch = nb_.initSearch(pbc, refsel);
157 dh.startFrame(frnr, fr.time);
158 for (size_t g = 0; g < sel_.size(); ++g)
160 const Selection &sel = pdata->parallelSelection(sel_[g]);
161 int nr = sel.posCount();
162 real frave = 0.0;
163 for (int i = 0; i < nr; ++i)
165 SelectionPosition p = sel.position(i);
166 frave += nbsearch.minimumDistance(p.x());
168 frave /= nr;
169 dh.setPoint(g, frave);
171 dh.finishFrame();
175 void
176 AnalysisTemplate::finishAnalysis(int /*nframes*/)
181 void
182 AnalysisTemplate::writeOutput()
184 // We print out the average of the mean distances for each group.
185 for (size_t g = 0; g < sel_.size(); ++g)
187 fprintf(stderr, "Average mean distance for '%s': %.3f nm\n",
188 sel_[g].name(), avem_->average(0, g));
192 /*! \brief
193 * The main function for the analysis template.
196 main(int argc, char *argv[])
198 return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<AnalysisTemplate>(argc, argv);