Be quieter when failing to detect GPUs
[gromacs.git] / share / template / template.cpp
blob5b82cd5fdd1ce5b4a156de7da49f269d31ebd55b
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(IOptionsContainer *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 : cutoff_(0.0)
79 registerAnalysisDataset(&data_, "avedist");
83 void
84 AnalysisTemplate::initOptions(IOptionsContainer *options,
85 TrajectoryAnalysisSettings *settings)
87 static const char *const desc[] = {
88 "This is a template for writing your own analysis tools for",
89 "GROMACS. The advantage of using GROMACS for this is that you",
90 "have access to all information in the topology, and your",
91 "program will be able to handle all types of coordinates and",
92 "trajectory files supported by GROMACS. In addition,",
93 "you get a lot of functionality for free from the trajectory",
94 "analysis library, including support for flexible dynamic",
95 "selections. Go ahead an try it![PAR]",
96 "To get started with implementing your own analysis program,",
97 "follow the instructions in the README file provided.",
98 "This template implements a simple analysis programs that calculates",
99 "average distances from a reference group to one or more",
100 "analysis groups."
103 settings->setHelpText(desc);
105 options->addOption(FileNameOption("o")
106 .filetype(eftPlot).outputFile()
107 .store(&fnDist_).defaultBasename("avedist")
108 .description("Average distances from reference group"));
110 options->addOption(SelectionOption("reference")
111 .store(&refsel_).required()
112 .description("Reference group to calculate distances from"));
113 options->addOption(SelectionOption("select")
114 .storeVector(&sel_).required().multiValue()
115 .description("Groups to calculate distances to"));
117 options->addOption(DoubleOption("cutoff").store(&cutoff_)
118 .description("Cutoff for distance calculation (0 = no cutoff)"));
120 settings->setFlag(TrajectoryAnalysisSettings::efRequireTop);
124 void
125 AnalysisTemplate::initAnalysis(const TrajectoryAnalysisSettings &settings,
126 const TopologyInformation & /*top*/)
128 nb_.setCutoff(cutoff_);
130 data_.setColumnCount(0, sel_.size());
132 avem_.reset(new AnalysisDataAverageModule());
133 data_.addModule(avem_);
135 if (!fnDist_.empty())
137 AnalysisDataPlotModulePointer plotm(
138 new AnalysisDataPlotModule(settings.plotSettings()));
139 plotm->setFileName(fnDist_);
140 plotm->setTitle("Average distance");
141 plotm->setXAxisIsTime();
142 plotm->setYLabel("Distance (nm)");
143 data_.addModule(plotm);
148 void
149 AnalysisTemplate::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
150 TrajectoryAnalysisModuleData *pdata)
152 AnalysisDataHandle dh = pdata->dataHandle(data_);
153 const Selection &refsel = pdata->parallelSelection(refsel_);
155 AnalysisNeighborhoodSearch nbsearch = nb_.initSearch(pbc, refsel);
156 dh.startFrame(frnr, fr.time);
157 for (size_t g = 0; g < sel_.size(); ++g)
159 const Selection &sel = pdata->parallelSelection(sel_[g]);
160 int nr = sel.posCount();
161 real frave = 0.0;
162 for (int i = 0; i < nr; ++i)
164 SelectionPosition p = sel.position(i);
165 frave += nbsearch.minimumDistance(p.x());
167 frave /= nr;
168 dh.setPoint(g, frave);
170 dh.finishFrame();
174 void
175 AnalysisTemplate::finishAnalysis(int /*nframes*/)
180 void
181 AnalysisTemplate::writeOutput()
183 // We print out the average of the mean distances for each group.
184 for (size_t g = 0; g < sel_.size(); ++g)
186 fprintf(stderr, "Average mean distance for '%s': %.3f nm\n",
187 sel_[g].name(), avem_->average(0, g));
191 /*! \brief
192 * The main function for the analysis template.
195 main(int argc, char *argv[])
197 return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<AnalysisTemplate>(argc, argv);