5c13d802c3843c21ee659dbc5f6dbc1451ffb0fc
[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 #include <gromacs/math/units.h>
41 #include <gromacs/topology/atomprop.h>
42
43
44 using namespace gmx;
45
46 /*! \brief
47  * Template class to serve as a basis for user analysis tools.
48  */
49 class SANS : public TrajectoryAnalysisModule
50 {
51     public:
52         SANS();
53
54         virtual void initOptions(IOptionsContainer          *options,
55                                  TrajectoryAnalysisSettings *settings);
56         virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
57                                   const TopologyInformation        &top);
58         virtual void initAfterFirstFrame(const TrajectoryAnalysisSettings &settings,
59                                          const t_trxframe                 &fr);
60
61
62         virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
63                                   TrajectoryAnalysisModuleData *pdata);
64
65         virtual void finishAnalysis(int nframes);
66         virtual void writeOutput();
67
68     private:
69         class ModuleData;
70
71         double Gaussian(double value, double sigma, RVec Rpos, RVec Rgrid, RVec GridSpacing);
72
73         double                           cutoff_;
74         double                           grid_;
75         double                           boxScale_;
76
77         IVec gridPoints_;
78         RVec gridSpacing_;
79
80         AnalysisNeighborhood             nb_;
81         const TopologyInformation       *top_;
82         std::vector < std::vector < std::vector<double>>> gausGrid_;
83         std::vector<double>               vdw_radius_;
84 };
85
86 SANS::SANS()
87     : cutoff_(1.0), grid_(0.05), boxScale_(1.0)
88 {
89 }
90
91 double SANS::Gaussian(double value, double sigma, RVec Rpos, RVec Rgrid, RVec GridSpacing)
92 {
93     //return 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*108.);
94     //return value*std::exp(-distance2(Rpos, Rgrid)/(2*sigma*sigma))/(2*M_PI*std::sqrt(2.0*M_PI)*sigma*sigma*sigma);
95     return value*std::exp(-distance2(Rpos, Rgrid)/(2*sigma*sigma))/(2*M_PI*std::sqrt(2.0*M_PI)*sigma*sigma*sigma)*GridSpacing[XX]*GridSpacing[YY]*GridSpacing[ZZ];
96     //return value*std::exp(-distance2(Rpos, Rgrid)/(sigma*sigma))/(2*M_PI*std::sqrt(2.0*M_PI)*sigma*sigma*sigma)*GridSpacing[XX]*GridSpacing[YY]*GridSpacing[ZZ];
97 }
98
99 void
100 SANS::initOptions(IOptionsContainer          *options,
101                   TrajectoryAnalysisSettings *settings)
102 {
103     static const char *const desc[] = {
104         "Analysis tool for finding molecular core."
105     };
106
107     // Add the descriptive text (program help text) to the options
108     settings->setHelpText(desc);
109     // Add option for cutoff constant
110     options->addOption(DoubleOption("cutoff")
111                            .store(&cutoff_)
112                            .description("cutoff for neighbour search"));
113     options->addOption(DoubleOption("grid")
114                            .store(&grid_)
115                            .description("grid spacing for gaus grid"));
116     options->addOption(DoubleOption("boxscale")
117                            .store(&boxScale_)
118                            .description("box scaling for gaus grid"));
119 }
120
121 void
122 SANS::initAnalysis(const TrajectoryAnalysisSettings  &settings,
123                    const TopologyInformation         &top)
124 {
125     nb_.setCutoff(cutoff_);
126     top_ = &top;
127     gmx_atomprop_t aps    = gmx_atomprop_init();
128     t_atoms       *atoms  = &(top.topology()->atoms);
129     real value;
130     for(int i=0; i<top.topology()->atoms.nr;i++) {
131         // Lookup the Van der Waals radius of this atom
132         int resnr = atoms->atom[i].resind;
133         if (gmx_atomprop_query(aps, epropVDW,
134                                *(atoms->resinfo[resnr].name),
135                                *(atoms->atomname[i]),
136                                &value))
137             {
138                 vdw_radius_.push_back(value);
139             }
140             else
141             {
142                 fprintf(stderr, "Could not determine VDW radius for %s-%s. Set to zero.\n",
143                                 *(atoms->resinfo[resnr].name),
144                                 *(atoms->atomname[i]));
145                     vdw_radius_.push_back(0.0);
146             }
147
148     }
149
150 }
151
152 void
153 SANS::initAfterFirstFrame(const TrajectoryAnalysisSettings &settings,
154                           const t_trxframe                 &fr)
155 {
156     // set gridpoints and grid spacing for orthogonal boxes only....
157     for (int i = 0; i < DIM; i++)
158     {
159         gridSpacing_[i] = grid_;
160         gridPoints_[i]  = static_cast<int>(std::ceil(fr.box[i][i]*boxScale_/grid_));
161     }
162     // prepare gausGrid
163     gausGrid_.resize(gridPoints_[XX], std::vector < std::vector < double>>(gridPoints_[YY], std::vector<double>(gridPoints_[ZZ])));
164 }
165
166 void
167 SANS::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
168                    TrajectoryAnalysisModuleData *pdata)
169 {
170     RVec GridSpacing(0.1, 0.1, 0.1);
171     AnalysisNeighborhoodSearch nbsearch = nb_.initSearch(pbc, AnalysisNeighborhoodPositions(fr.x, fr.natoms));
172
173     t_topology                *top = top_->topology();
174
175     fprintf(stderr, "Box dimensions:\n");
176     for (int i = 0; i < DIM; i++)
177     {
178         fprintf(stderr, "[ ");
179         for (int j = 0; j < DIM; j++)
180         {
181             fprintf(stderr, " %f ", fr.box[i][j]);
182         }
183         fprintf(stderr, " ]\n");
184     }
185     // main loop over grid points
186     for (int i = 0; i < gridPoints_[XX]; i++)
187     {
188         for (int j = 0; j < gridPoints_[YY]; j++)
189         {
190             for (int k = 0; k < gridPoints_[ZZ]; k++)
191             {
192                 RVec                           point(i*gridSpacing_[XX], j*gridSpacing_[YY], k*gridSpacing_[ZZ]);
193                 AnalysisNeighborhoodPairSearch pairSearch = nbsearch.startPairSearch(point.as_vec());
194                 AnalysisNeighborhoodPair       pair;
195                 while (pairSearch.findNextPair(&pair))
196                 {
197                     //fprintf(stderr,"point = [ %f %f %f ]\n", point[XX],point[YY],point[ZZ]);
198                     // fprintf(stderr,"Index %d\n", pair.refIndex());
199                     // fprintf(stderr,"dx =  (%f, %f, %f)\n", pair.dx()[XX], pair.dx()[YY], pair.dx()[ZZ] );
200                     // fprintf(stderr,"ref =  (%f, %f, %f)\n", fr.x[pair.refIndex()][XX], fr.x[pair.refIndex()][YY], fr.x[pair.refIndex()][ZZ]);
201                     // fprintf(stderr,"ref_dx =  (%f, %f, %f)\n", GridSpacing[XX]+pair.dx()[XX], GridSpacing[YY]+pair.dx()[YY], GridSpacing[ZZ]+pair.dx()[ZZ] );
202                     RVec refPos(point[XX]+pair.dx()[XX], point[YY]+pair.dx()[YY], point[ZZ]+pair.dx()[ZZ]);
203                     //fprintf(stderr,"refPos = [ %f %f %f ]\n", refPos[XX], refPos[YY], refPos[ZZ]);
204                     //fprintf(stderr, "Mass = %f, Radius = %f\n", top->atoms.atom[pair.refIndex()].m, top->atomtypes.radius[top->atoms.atom[pair.refIndex()].type]);
205                     //fprintf(stderr, "Distance^2 = %3.8f\n", pair.distance2());
206                     //fprintf(stderr, "Gausian = %3.8f\n",Gaussian(top->atoms.atom[pair.refIndex()].m, 0.2, refPos, point, gridSpacing_)*AMU/((NANO*NANO*NANO)*(gridSpacing_[XX]*gridSpacing_[YY]*gridSpacing_[ZZ])));
207                     //gausGrid_[i][j][k] += Gaussian(top->atoms.atom[pair.refIndex()].m, 0.1, refPos, point, gridSpacing_)*AMU/((NANO*NANO*NANO)*(gridSpacing_[XX]*gridSpacing_[YY]*gridSpacing_[ZZ]));
208                     gausGrid_[i][j][k] += Gaussian(top->atoms.atom[pair.refIndex()].m, vdw_radius_[pair.refIndex()], refPos, point, gridSpacing_)*AMU/((NANO*NANO*NANO)*(gridSpacing_[XX]*gridSpacing_[YY]*gridSpacing_[ZZ]));
209                 }
210                 fprintf(stderr, "GausGrid[%d, %d, %d] = %f\n", i, j, k, gausGrid_[i][j][k]);
211             }
212         }
213     }
214     // only for orthogonal boxes...
215     fprintf(stderr, "NX = %d, NY = %d, NZ = %d\n", gridPoints_[XX], gridPoints_[YY], gridPoints_[ZZ]);
216 }
217
218 void
219 SANS::finishAnalysis(int nframes)
220 {
221 }
222
223 void
224 SANS::writeOutput()
225 {
226     FILE *fo;
227     // Construct opendx grid
228     fo = fopen("data.dx","w");
229     //object 1 class gridpositions counts  140 140 140
230     //origin 0.000000 0.000000 0.000000
231     //delta  0.500000 0.000000 0.000000
232     //delta  0.000000 0.500000 0.000000
233     //delta  0.000000 0.000000 0.500000
234     //object 2 class gridconnections counts  140 140 140
235     //object 3 class array type "double" rank 0 items 2744000 data follows
236     fprintf(fo, "object 1 class gridpositions counts %d %d %d\n", gridPoints_[XX],gridPoints_[YY],gridPoints_[ZZ]);
237     fprintf(fo, "origin 0.000000 0.000000 0.000000\n");
238     fprintf(fo, "delta %3.6f %3.6f %3.6f\n", gridSpacing_[XX]*NM2A, 0.0, 0.0);
239     fprintf(fo, "delta %3.6f %3.6f %3.6f\n", 0.0, gridSpacing_[YY]*NM2A, 0.0);
240     fprintf(fo, "delta %3.6f %3.6f %3.6f\n", 0.0, 0.0, gridSpacing_[ZZ]*NM2A);
241     fprintf(fo, "object 2 class gridconnections counts  %d %d %d\n", gridPoints_[XX],gridPoints_[YY],gridPoints_[ZZ]);
242     fprintf(fo, "object 3 class array type \"double\" rank 0 items %d data follows\n", gridPoints_[XX]*gridPoints_[YY]*gridPoints_[ZZ]);
243     // now dump data in ordered mode
244     //    u(0,0,0) u(0,0,1) u(0,0,2)
245     //    ...
246     //    u(0,0,nz-3) u(0,0,nz-2) u(0,0,nz-1)
247     //    u(0,1,0) u(0,1,1) u(0,1,2)
248     //    ...
249     //    u(0,1,nz-3) u(0,1,nz-2) u(0,1,nz-1)
250     //    ...
251     //    u(0,ny-1,nz-3) u(0,ny-1,nz-2) u(0,ny-1,nz-1)
252     //    u(1,0,0) u(1,0,1) u(1,0,2)
253     //    ...
254     int NR = 3;
255     int n = 0;
256     for (int i = 0; i < gridPoints_[XX]; i++)
257     {
258         for (int j = 0; j < gridPoints_[YY]; j++)
259         {
260             for (int k = 0; k < gridPoints_[ZZ]; k++)
261             {
262                 fprintf(fo, "%3.8f ", gausGrid_[i][j][k]);
263                 n++;
264                 if (n == NR) {
265                     fprintf(fo, "\n");
266                     n = 0;
267                 }
268             }
269         }
270     }
271     if (n != NR) {
272         fprintf(fo, "\n");
273     }
274     //attribute "dep" string "positions"
275     //object "density" class field
276     //component "positions" value 1
277     //component "connections" value 2
278     //component "data" value 3
279     fprintf(fo, "attribute \"dep\" string \"positions\"\n");
280     fprintf(fo, "object \"density\" class field\n");
281     fprintf(fo, "component \"positions\" value 1\n");
282     fprintf(fo, "component \"connections\" value 2\n");
283     fprintf(fo, "component \"data\" value 3\n");
284     fclose(fo);
285 }
286
287 /*! \brief
288  * The main function for the analysis template.
289  */
290 int
291 main(int argc, char *argv[])
292 {
293     return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<SANS>(argc, argv);
294 }