Merge branch 'master' into pygromacs
[alexxy/gromacs.git] / src / gromacs / trajectoryanalysis / cmdlinerunner.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,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.
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 gmx::TrajectoryAnalysisCommandLineRunner.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_trajectoryanalysis
41  */
42 #include "gmxpre.h"
43
44 #include "cmdlinerunner.h"
45
46 #include "gromacs/analysisdata/paralleloptions.h"
47 #include "gromacs/commandline/cmdlinehelpcontext.h"
48 #include "gromacs/commandline/cmdlinehelpwriter.h"
49 #include "gromacs/commandline/cmdlinemodule.h"
50 #include "gromacs/commandline/cmdlinemodulemanager.h"
51 #include "gromacs/commandline/cmdlineparser.h"
52 #include "gromacs/fileio/trx.h"
53 #include "gromacs/options/filenameoptionmanager.h"
54 #include "gromacs/options/options.h"
55 #include "gromacs/pbcutil/pbc.h"
56 #include "gromacs/selection/selectioncollection.h"
57 #include "gromacs/selection/selectionoptionmanager.h"
58 #include "gromacs/trajectoryanalysis/analysismodule.h"
59 #include "gromacs/trajectoryanalysis/analysissettings.h"
60 #include "gromacs/utility/exceptions.h"
61 #include "gromacs/utility/filestream.h"
62 #include "gromacs/utility/gmxassert.h"
63
64 #include "runnercommon.h"
65
66 namespace gmx
67 {
68
69 /********************************************************************
70  * TrajectoryAnalysisCommandLineRunner::Impl
71  */
72
73 class TrajectoryAnalysisCommandLineRunner::Impl
74 {
75     public:
76         class RunnerCommandLineModule;
77
78         Impl(TrajectoryAnalysisModule *module);
79         ~Impl();
80
81         void parseOptions(TrajectoryAnalysisSettings *settings,
82                           TrajectoryAnalysisRunnerCommon *common,
83                           SelectionCollection *selections,
84                           int *argc, char *argv[], bool full=true);
85
86         TrajectoryAnalysisModule *module_;
87         bool                      bUseDefaultGroups_;
88         int                       debugLevel_;
89 };
90
91
92 TrajectoryAnalysisCommandLineRunner::Impl::Impl(
93         TrajectoryAnalysisModule *module)
94     : module_(module), bUseDefaultGroups_(true), debugLevel_(0)
95 {
96 }
97
98
99 TrajectoryAnalysisCommandLineRunner::Impl::~Impl()
100 {
101 }
102
103
104 void
105 TrajectoryAnalysisCommandLineRunner::Impl::parseOptions(
106         TrajectoryAnalysisSettings *settings,
107         TrajectoryAnalysisRunnerCommon *common,
108         SelectionCollection *selections,
109         int *argc, char *argv[], bool full)
110 {
111     FileNameOptionManager  fileoptManager;
112     SelectionOptionManager seloptManager(selections);
113     Options                options(NULL, NULL);
114     Options                moduleOptions(module_->name(), module_->description());
115     Options                commonOptions("common", "Common analysis control");
116     Options                selectionOptions("selection", "Common selection control");
117
118     options.addManager(&fileoptManager);
119     options.addManager(&seloptManager);
120     if (full) {
121         options.addSubSection(&commonOptions);
122     }
123     options.addSubSection(&selectionOptions);
124     options.addSubSection(&moduleOptions);
125
126     module_->initOptions(&moduleOptions, settings);
127     if (full) {
128         common->initOptions(&commonOptions);
129     }
130     selections->initOptions(&selectionOptions);
131
132     {
133         CommandLineParser  parser(&options);
134         // TODO: Print the help if user provides an invalid option?
135         // Or just add a message advising the user to invoke the help?
136         parser.parse(argc, argv);
137         common->scaleTimeOptions(&options);
138         options.finish();
139     }
140
141     if (full) {
142         common->optionsFinished(&commonOptions);
143     }
144     module_->optionsFinished(&moduleOptions, settings);
145
146     common->initIndexGroups(selections, bUseDefaultGroups_);
147
148     const bool bInteractive = StandardInputStream::instance().isInteractive();
149     seloptManager.parseRequestedFromStdin(bInteractive);
150
151     common->doneIndexGroups(selections);
152     common->initTopology(selections);
153
154     selections->compile();
155 }
156
157
158 /********************************************************************
159  * TrajectoryAnalysisCommandLineRunner
160  */
161
162 TrajectoryAnalysisCommandLineRunner::TrajectoryAnalysisCommandLineRunner(
163         TrajectoryAnalysisModule *module)
164     : impl_(new Impl(module))
165 {
166 }
167
168
169 TrajectoryAnalysisCommandLineRunner::~TrajectoryAnalysisCommandLineRunner()
170 {
171 }
172
173
174 void
175 TrajectoryAnalysisCommandLineRunner::setUseDefaultGroups(bool bUseDefaults)
176 {
177     impl_->bUseDefaultGroups_ = bUseDefaults;
178 }
179
180
181 void
182 TrajectoryAnalysisCommandLineRunner::setSelectionDebugLevel(int debuglevel)
183 {
184     impl_->debugLevel_ = debuglevel;
185 }
186
187
188 int
189 TrajectoryAnalysisCommandLineRunner::run(int argc, char *argv[])
190 {
191     TrajectoryAnalysisModule *module = impl_->module_;
192
193     SelectionCollection       selections;
194     selections.setDebugLevel(impl_->debugLevel_);
195
196     TrajectoryAnalysisSettings      settings;
197     TrajectoryAnalysisRunnerCommon  common(&settings);
198
199     impl_->parseOptions(&settings, &common, &selections, &argc, argv);
200
201     const TopologyInformation &topology = common.topologyInformation();
202     module->initAnalysis(settings, topology);
203
204     TrajectoryAnalysisModule::Batch batch = module->getBatch();
205     std::vector<SelectionCollection*> batchSelections;
206     std::vector<Impl*> impls;
207     for (size_t i = 0; i < batch.size(); i++) {
208         TrajectoryAnalysisModule *bmodule = batch[i];
209         batchSelections.push_back(new SelectionCollection());
210         impls.push_back(new Impl(bmodule));
211         std::vector<char*> modArgv = module->getArgv(i);
212         int modArgc = modArgv.size();
213         impls.back()->parseOptions(&settings, &common, batchSelections.back(), &modArgc, modArgv.data(), false);
214
215         batch[i]->initAnalysis(settings, topology);
216     }
217
218     // Load first frame.
219     common.initFirstFrame();
220     module->initAfterFirstFrame(settings, common.frame());
221
222     t_pbc  pbc;
223     t_pbc *ppbc = settings.hasPBC() ? &pbc : NULL;
224
225     int    nframes = 0;
226     AnalysisDataParallelOptions         dataOptions;
227     TrajectoryAnalysisModuleDataPointer pdata(
228             module->startFrames(dataOptions, selections));
229
230     std::vector<AnalysisDataParallelOptions> batchOptions;
231     std::vector<TrajectoryAnalysisModuleDataPointer> batchDataPointers;
232     for (size_t i = 0; i < batch.size(); i++) {
233         batch[i]->initAfterFirstFrame(settings, common.frame());
234
235         batchOptions.push_back(AnalysisDataParallelOptions());
236         batchDataPointers.push_back(batch[i]->startFrames(
237             batchOptions.back(), *batchSelections[i]));
238     }
239
240     do
241     {
242         common.initFrame();
243         t_trxframe &frame = common.frame();
244         if (ppbc != NULL)
245         {
246             set_pbc(ppbc, topology.ePBC(), frame.box);
247         }
248
249         for (size_t i = 0; i < batch.size(); i++) {
250             batchSelections[i]->evaluate(&frame, ppbc);
251             batch[i]->analyzeFrame(nframes, frame, ppbc, batchDataPointers[i].get());
252             batch[i]->finishFrameSerial(nframes);
253         }
254         selections.evaluate(&frame, ppbc);
255         module->analyzeFrame(nframes, frame, ppbc, pdata.get());
256         module->finishFrameSerial(nframes);
257
258         ++nframes;
259     }
260     while (common.readNextFrame());
261     for (size_t i = 0; i < batch.size(); i++) {
262         batch[i]->finishFrames(batchDataPointers[i].get());
263         if (batchDataPointers[i].get() != NULL) {
264             batchDataPointers[i]->finish();
265         }
266         batchDataPointers[i].reset();
267     }
268
269     module->finishFrames(pdata.get());
270     if (pdata.get() != NULL)
271     {
272         pdata->finish();
273     }
274     pdata.reset();
275
276     if (common.hasTrajectory())
277     {
278         fprintf(stderr, "Analyzed %d frames, last time %.3f\n",
279                 nframes, common.frame().time);
280     }
281     else
282     {
283         fprintf(stderr, "Analyzed topology coordinates\n");
284     }
285
286     // Restore the maximal groups for dynamic selections.
287     for (size_t i = 0; i < batch.size(); i++) {
288         batchSelections[i]->evaluateFinal(nframes);
289         batch[i]->finishAnalysis(nframes);
290         batch[i]->writeOutput();
291
292         delete batchSelections[i];
293         delete impls[i];
294     }
295
296     selections.evaluateFinal(nframes);
297
298     module->finishAnalysis(nframes);
299     module->writeOutput();
300
301     return 0;
302 }
303
304
305 void
306 TrajectoryAnalysisCommandLineRunner::writeHelp(const CommandLineHelpContext &context)
307 {
308     // TODO: This method duplicates some code from run().
309     // See how to best refactor it to share the common code.
310     SelectionCollection             selections;
311     TrajectoryAnalysisSettings      settings;
312     TrajectoryAnalysisRunnerCommon  common(&settings);
313
314     SelectionOptionManager          seloptManager(&selections);
315     Options                         options(NULL, NULL);
316     Options                         moduleOptions(impl_->module_->name(), impl_->module_->description());
317     Options                         commonOptions("common", "Common analysis control");
318     Options                         selectionOptions("selection", "Common selection control");
319
320     options.addManager(&seloptManager);
321     options.addSubSection(&commonOptions);
322     options.addSubSection(&selectionOptions);
323     options.addSubSection(&moduleOptions);
324
325     impl_->module_->initOptions(&moduleOptions, &settings);
326     common.initOptions(&commonOptions);
327     selections.initOptions(&selectionOptions);
328
329     CommandLineHelpWriter(options)
330         .setHelpText(settings.helpText())
331         .setTimeUnitString(settings.timeUnitManager().timeUnitAsString())
332         .writeHelp(context);
333 }
334
335
336 /*! \internal \brief
337  * Command line module for a trajectory analysis module.
338  *
339  * \ingroup module_trajectoryanalysis
340  */
341 class TrajectoryAnalysisCommandLineRunner::Impl::RunnerCommandLineModule
342     : public ICommandLineModule
343 {
344     public:
345         /*! \brief
346          * Constructs a module.
347          *
348          * \param[in] name         Name for the module.
349          * \param[in] description  One-line description for the module.
350          * \param[in] factory      Factory method to create the analysis module.
351          *
352          * Does not throw.  This is important for correct implementation of
353          * runAsMain().
354          */
355         RunnerCommandLineModule(const char *name, const char *description,
356                                 ModuleFactoryMethod factory)
357             : name_(name), description_(description), hasFunction_(true), factory_(factory), functor_(NULL)
358         {
359         }
360
361         /*! \brief
362          * Overloaded constructor accepting a functor instead of function pointer.
363          */
364         RunnerCommandLineModule(const char *name, const char *description,
365                                 ModuleFactoryFunctor *factory)
366             : name_(name), description_(description), hasFunction_(false), factory_(NULL), functor_(factory)
367         {
368         }
369
370         virtual const char *name() const { return name_; }
371         virtual const char *shortDescription() const { return description_; };
372
373         virtual void init(CommandLineModuleSettings *settings);
374         virtual int run(int argc, char *argv[]);
375         virtual void writeHelp(const CommandLineHelpContext &context) const;
376
377     private:
378         const char            *name_;
379         const char            *description_;
380         bool                   hasFunction_;
381         ModuleFactoryMethod    factory_;
382         ModuleFactoryFunctor   *functor_;
383
384         GMX_DISALLOW_COPY_AND_ASSIGN(RunnerCommandLineModule);
385 };
386
387 void TrajectoryAnalysisCommandLineRunner::Impl::RunnerCommandLineModule::init(
388         CommandLineModuleSettings * /*settings*/)
389 {
390 }
391
392 int TrajectoryAnalysisCommandLineRunner::Impl::RunnerCommandLineModule::run(
393         int argc, char *argv[])
394 {
395     TrajectoryAnalysisModulePointer     module(hasFunction_? factory_() : (*functor_)());
396     TrajectoryAnalysisCommandLineRunner runner(module.get());
397     return runner.run(argc, argv);
398 }
399
400 void TrajectoryAnalysisCommandLineRunner::Impl::RunnerCommandLineModule::writeHelp(
401         const CommandLineHelpContext &context) const
402 {
403     TrajectoryAnalysisModulePointer     module(hasFunction_? factory_() : (*functor_)());
404     TrajectoryAnalysisCommandLineRunner runner(module.get());
405     runner.writeHelp(context);
406 }
407
408 // static
409 int
410 TrajectoryAnalysisCommandLineRunner::runAsMain(
411         int argc, char *argv[], ModuleFactoryMethod factory)
412 {
413     Impl::RunnerCommandLineModule module(NULL, NULL, factory);
414     return CommandLineModuleManager::runAsMainSingleModule(argc, argv, &module);
415 }
416
417 // static
418 int
419 TrajectoryAnalysisCommandLineRunner::runAsMain(
420         int argc, char *argv[], ModuleFactoryFunctor *factory)
421 {
422     Impl::RunnerCommandLineModule module(NULL, NULL, factory);
423     return CommandLineModuleManager::runAsMainSingleModule(argc, argv, &module);
424 }
425
426 // static
427 void
428 TrajectoryAnalysisCommandLineRunner::registerModule(
429         CommandLineModuleManager *manager, const char *name,
430         const char *description, ModuleFactoryMethod factory)
431 {
432     CommandLineModulePointer module(
433             new Impl::RunnerCommandLineModule(name, description, factory));
434     manager->addModule(move(module));
435 }
436
437 } // namespace gmx