bfd4d0efd95f8477239eb348d3ca6aaa49534dfb
[alexxy/gromacs-sans.git] / src / sans.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 #include <iostream>
36
37 #include <gromacs/trajectoryanalysis.h>
38 #include <gromacs/selection/nbsearch.h>
39 #include <gromacs/math/vec.h>
40
41 using namespace gmx;
42
43 /*! \brief
44  * Template class to serve as a basis for user analysis tools.
45  */
46 class SANS : public TrajectoryAnalysisModule
47 {
48     public:
49         SANS();
50
51         virtual void initOptions(IOptionsContainer          *options,
52                                  TrajectoryAnalysisSettings *settings);
53         virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
54                                   const TopologyInformation        &top);
55
56         virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
57                                   TrajectoryAnalysisModuleData *pdata);
58
59         virtual void finishAnalysis(int nframes);
60         virtual void writeOutput();
61
62     private:
63         class ModuleData;
64
65         double Gaussian(double value, double sigma, RVec Rpos, RVec Rgrid, RVec GridSpacing);
66
67         std::string                      fnNdx_;
68         double                           cutoff_;
69         double                           grid_;
70
71         AnalysisNeighborhood             nb_;
72         const TopologyInformation          *top_;
73
74 };
75
76 SANS::SANS()
77     : cutoff_(1.0), grid_(0.05)
78 {
79 }
80
81 double SANS::Gaussian(double value, double sigma, RVec Rpos, RVec Rgrid, RVec GridSpacing)
82 {
83         return 2.*value*M_PI*M_PI*GridSpacing[XX]*GridSpacing[YY]*GridSpacing[ZZ]*exp(-distance2(Rpos, Rgrid)/(2*sigma*sigma))/(std::sqrt(2.0)*sigma*sigma*sigma*216.);
84 }
85
86 void
87 SANS::initOptions(IOptionsContainer          *options,
88                               TrajectoryAnalysisSettings *settings)
89 {
90     static const char *const desc[] = {
91         "Analysis tool for finding molecular core."
92     };
93
94     // Add the descriptive text (program help text) to the options
95     settings->setHelpText(desc);
96 }
97
98 void
99 SANS::initAnalysis(const TrajectoryAnalysisSettings &settings,
100                                const TopologyInformation         & top)
101 {
102     nb_.setCutoff(cutoff_);
103     top_ = &top;
104 }
105
106 void
107 SANS::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
108                                TrajectoryAnalysisModuleData *pdata)
109 {
110     RVec GridSpacing(0.1, 0.1, 0.1);
111     AnalysisNeighborhoodSearch nbsearch = nb_.initSearch(pbc, AnalysisNeighborhoodPositions(fr.x, fr.natoms));
112     AnalysisNeighborhoodPairSearch pairSearch = nbsearch.startPairSearch(GridSpacing.as_vec());
113     AnalysisNeighborhoodPair pair;
114     t_topology *top = top_->topology();
115     while (pairSearch.findNextPair(&pair))
116       {
117           fprintf(stderr,"Index %d\n", pair.refIndex());
118           fprintf(stderr,"dx =  (%f, %f, %f)\n", pair.dx()[XX], pair.dx()[YY], pair.dx()[ZZ] );
119           fprintf(stderr,"ref =  (%f, %f, %f)\n", fr.x[pair.refIndex()][XX], fr.x[pair.refIndex()][YY], fr.x[pair.refIndex()][ZZ]);
120           fprintf(stderr,"ref_dx =  (%f, %f, %f)\n", GridSpacing[XX]+pair.dx()[XX], GridSpacing[YY]+pair.dx()[YY], GridSpacing[ZZ]+pair.dx()[ZZ] );
121
122       }
123 }
124
125 void
126 SANS::finishAnalysis(int nframes)
127 {
128 }
129
130 void
131 SANS::writeOutput()
132 {
133
134 }
135
136 /*! \brief
137  * The main function for the analysis template.
138  */
139 int
140 main(int argc, char *argv[])
141 {
142     return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<SANS>(argc, argv);
143 }