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