6838d272a3fa2cb15d8d11fda3578300b31d0b3c
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_sans.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, 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
36 #include "gmxpre.h"
37
38 #include "config.h"
39
40 #include "gromacs/legacyheaders/typedefs.h"
41 #include "gromacs/legacyheaders/macros.h"
42 #include "gromacs/math/vec.h"
43 #include "gromacs/legacyheaders/copyrite.h"
44 #include "gromacs/topology/index.h"
45 #include "gstat.h"
46 #include "gmx_ana.h"
47 #include "nsfactor.h"
48
49 #include "gromacs/commandline/pargs.h"
50 #include "gromacs/fileio/tpxio.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "gromacs/fileio/xvgr.h"
53 #include "gromacs/pbcutil/rmpbc.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/utility/gmxomp.h"
57 #include "gromacs/utility/smalloc.h"
58
59 int gmx_sans(int argc, char *argv[])
60 {
61     const char          *desc[] = {
62         "[THISMODULE] computes SANS spectra using Debye formula.",
63         "It currently uses topology file (since it need to assigne element for each atom).",
64         "[PAR]",
65         "Parameters:[PAR]"
66         "[TT]-pr[tt] Computes normalized g(r) function averaged over trajectory[PAR]",
67         "[TT]-prframe[tt] Computes normalized g(r) function for each frame[PAR]",
68         "[TT]-sq[tt] Computes SANS intensity curve averaged over trajectory[PAR]",
69         "[TT]-sqframe[tt] Computes SANS intensity curve for each frame[PAR]",
70         "[TT]-startq[tt] Starting q value in nm[PAR]",
71         "[TT]-endq[tt] Ending q value in nm[PAR]",
72         "[TT]-qstep[tt] Stepping in q space[PAR]",
73         "Note: When using Debye direct method computational cost increases as",
74         "1/2 * N * (N - 1) where N is atom number in group of interest.",
75         "[PAR]",
76         "WARNING: If sq or pr specified this tool can produce large number of files! Up to two times larger than number of frames!"
77     };
78     static gmx_bool      bPBC     = TRUE;
79     static gmx_bool      bNORM    = FALSE;
80     static real          binwidth = 0.2, grid = 0.05; /* bins shouldnt be smaller then smallest bond (~0.1nm) length */
81     static real          start_q  = 0.0, end_q = 2.0, q_step = 0.01;
82     static real          mcover   = -1;
83     static unsigned int  seed     = 0;
84     static int           nthreads = -1;
85
86     static const char   *emode[]   = { NULL, "direct", "mc", NULL };
87     static const char   *emethod[] = { NULL, "debye", "fft", NULL };
88
89     gmx_neutron_atomic_structurefactors_t    *gnsf;
90     gmx_sans_t                               *gsans;
91
92 #define NPA asize(pa)
93
94     t_pargs                               pa[] = {
95         { "-bin", FALSE, etREAL, {&binwidth},
96           "[HIDDEN]Binwidth (nm)" },
97         { "-mode", FALSE, etENUM, {emode},
98           "Mode for sans spectra calculation" },
99         { "-mcover", FALSE, etREAL, {&mcover},
100           "Monte-Carlo coverage should be -1(default) or (0,1]"},
101         { "-method", FALSE, etENUM, {emethod},
102           "[HIDDEN]Method for sans spectra calculation" },
103         { "-pbc", FALSE, etBOOL, {&bPBC},
104           "Use periodic boundary conditions for computing distances" },
105         { "-grid", FALSE, etREAL, {&grid},
106           "[HIDDEN]Grid spacing (in nm) for FFTs" },
107         {"-startq", FALSE, etREAL, {&start_q},
108          "Starting q (1/nm) "},
109         {"-endq", FALSE, etREAL, {&end_q},
110          "Ending q (1/nm)"},
111         { "-qstep", FALSE, etREAL, {&q_step},
112           "Stepping in q (1/nm)"},
113         { "-seed",     FALSE, etINT,  {&seed},
114           "Random seed for Monte-Carlo"},
115 #ifdef GMX_OPENMP
116         { "-nt",  FALSE, etINT, {&nthreads},
117           "Number of threads to start"},
118 #endif
119     };
120     FILE                                 *fp;
121     const char                           *fnTPX, *fnNDX, *fnTRX, *fnDAT = NULL;
122     t_trxstatus                          *status;
123     t_topology                           *top  = NULL;
124     t_atom                               *atom = NULL;
125     gmx_rmpbc_t                           gpbc = NULL;
126     gmx_bool                              bTPX;
127     gmx_bool                              bFFT = FALSE, bDEBYE = FALSE;
128     gmx_bool                              bMC  = FALSE;
129     int                                   ePBC = -1;
130     matrix                                box;
131     char                                  title[STRLEN];
132     rvec                                 *x;
133     int                                   natoms;
134     real                                  t;
135     char                                **grpname = NULL;
136     atom_id                              *index   = NULL;
137     int                                   isize;
138     int                                   i, j;
139     char                                 *hdr            = NULL;
140     char                                 *suffix         = NULL;
141     t_filenm                             *fnmdup         = NULL;
142     gmx_radial_distribution_histogram_t  *prframecurrent = NULL, *pr = NULL;
143     gmx_static_structurefactor_t         *sqframecurrent = NULL, *sq = NULL;
144     output_env_t                          oenv;
145
146 #define NFILE asize(fnm)
147
148     t_filenm   fnm[] = {
149         { efTPX,  "-s",       NULL,       ffREAD },
150         { efTRX,  "-f",       NULL,       ffREAD },
151         { efNDX,  NULL,       NULL,       ffOPTRD },
152         { efDAT,  "-d",       "nsfactor", ffOPTRD },
153         { efXVG,  "-pr",      "pr",       ffWRITE },
154         { efXVG,  "-sq",       "sq",      ffWRITE },
155         { efXVG,  "-prframe", "prframe",  ffOPTWR },
156         { efXVG,  "-sqframe", "sqframe",  ffOPTWR }
157     };
158
159     nthreads = gmx_omp_get_max_threads();
160
161     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT,
162                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
163     {
164         return 0;
165     }
166
167     /* check that binwidth not smaller than smallers distance */
168     check_binwidth(binwidth);
169     check_mcover(mcover);
170
171     /* setting number of omp threads globaly */
172     gmx_omp_set_num_threads(nthreads);
173
174     /* Now try to parse opts for modes */
175     switch (emethod[0][0])
176     {
177         case 'd':
178             bDEBYE = TRUE;
179             switch (emode[0][0])
180             {
181                 case 'd':
182                     bMC = FALSE;
183                     break;
184                 case 'm':
185                     bMC = TRUE;
186                     break;
187                 default:
188                     break;
189             }
190             break;
191         case 'f':
192             bFFT = TRUE;
193             break;
194         default:
195             break;
196     }
197
198     if (bDEBYE)
199     {
200         if (bMC)
201         {
202             fprintf(stderr, "Using Monte Carlo Debye method to calculate spectrum\n");
203         }
204         else
205         {
206             fprintf(stderr, "Using direct Debye method to calculate spectrum\n");
207         }
208     }
209     else if (bFFT)
210     {
211         gmx_fatal(FARGS, "FFT method not implemented!");
212     }
213     else
214     {
215         gmx_fatal(FARGS, "Unknown combination for mode and method!");
216     }
217
218     /* Try to read files */
219     fnDAT = ftp2fn(efDAT, NFILE, fnm);
220     fnTPX = ftp2fn(efTPX, NFILE, fnm);
221     fnTRX = ftp2fn(efTRX, NFILE, fnm);
222
223     gnsf = gmx_neutronstructurefactors_init(fnDAT);
224     fprintf(stderr, "Read %d atom names from %s with neutron scattering parameters\n\n", gnsf->nratoms, fnDAT);
225
226     snew(top, 1);
227     snew(grpname, 1);
228     snew(index, 1);
229
230     bTPX = read_tps_conf(fnTPX, title, top, &ePBC, &x, NULL, box, TRUE);
231
232     printf("\nPlease select group for SANS spectra calculation:\n");
233     get_index(&(top->atoms), ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, grpname);
234
235     gsans = gmx_sans_init(top, gnsf);
236
237     /* Prepare reference frame */
238     if (bPBC)
239     {
240         gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
241         gmx_rmpbc(gpbc, top->atoms.nr, box, x);
242     }
243
244     natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box);
245     if (natoms != top->atoms.nr)
246     {
247         fprintf(stderr, "\nWARNING: number of atoms in tpx (%d) and trajectory (%d) do not match\n", natoms, top->atoms.nr);
248     }
249
250     do
251     {
252         if (bPBC)
253         {
254             gmx_rmpbc(gpbc, top->atoms.nr, box, x);
255         }
256         /* allocate memory for pr */
257         if (pr == NULL)
258         {
259             /* in case its first frame to read */
260             snew(pr, 1);
261         }
262         /*  realy calc p(r) */
263         prframecurrent = calc_radial_distribution_histogram(gsans, x, box, index, isize, binwidth, bMC, bNORM, mcover, seed);
264         /* copy prframecurrent -> pr and summ up pr->gr[i] */
265         /* allocate and/or resize memory for pr->gr[i] and pr->r[i] */
266         if (pr->gr == NULL)
267         {
268             /* check if we use pr->gr first time */
269             snew(pr->gr, prframecurrent->grn);
270             snew(pr->r, prframecurrent->grn);
271         }
272         else
273         {
274             /* resize pr->gr and pr->r if needed to preven overruns */
275             if (prframecurrent->grn > pr->grn)
276             {
277                 srenew(pr->gr, prframecurrent->grn);
278                 srenew(pr->r, prframecurrent->grn);
279             }
280         }
281         pr->grn      = prframecurrent->grn;
282         pr->binwidth = prframecurrent->binwidth;
283         /* summ up gr and fill r */
284         for (i = 0; i < prframecurrent->grn; i++)
285         {
286             pr->gr[i] += prframecurrent->gr[i];
287             pr->r[i]   = prframecurrent->r[i];
288         }
289         /* normalize histo */
290         normalize_probability(prframecurrent->grn, prframecurrent->gr);
291         /* convert p(r) to sq */
292         sqframecurrent = convert_histogram_to_intensity_curve(prframecurrent, start_q, end_q, q_step);
293         /* print frame data if needed */
294         if (opt2fn_null("-prframe", NFILE, fnm))
295         {
296             snew(hdr, 25);
297             snew(suffix, GMX_PATH_MAX);
298             /* prepare header */
299             sprintf(hdr, "g(r), t = %f", t);
300             /* prepare output filename */
301             fnmdup = dup_tfn(NFILE, fnm);
302             sprintf(suffix, "-t%.2f", t);
303             add_suffix_to_output_names(fnmdup, NFILE, suffix);
304             fp = xvgropen(opt2fn_null("-prframe", NFILE, fnmdup), hdr, "Distance (nm)", "Probability", oenv);
305             for (i = 0; i < prframecurrent->grn; i++)
306             {
307                 fprintf(fp, "%10.6f%10.6f\n", prframecurrent->r[i], prframecurrent->gr[i]);
308             }
309             done_filenms(NFILE, fnmdup);
310             fclose(fp);
311             sfree(hdr);
312             sfree(suffix);
313             sfree(fnmdup);
314         }
315         if (opt2fn_null("-sqframe", NFILE, fnm))
316         {
317             snew(hdr, 25);
318             snew(suffix, GMX_PATH_MAX);
319             /* prepare header */
320             sprintf(hdr, "I(q), t = %f", t);
321             /* prepare output filename */
322             fnmdup = dup_tfn(NFILE, fnm);
323             sprintf(suffix, "-t%.2f", t);
324             add_suffix_to_output_names(fnmdup, NFILE, suffix);
325             fp = xvgropen(opt2fn_null("-sqframe", NFILE, fnmdup), hdr, "q (nm^-1)", "s(q)/s(0)", oenv);
326             for (i = 0; i < sqframecurrent->qn; i++)
327             {
328                 fprintf(fp, "%10.6f%10.6f\n", sqframecurrent->q[i], sqframecurrent->s[i]);
329             }
330             done_filenms(NFILE, fnmdup);
331             fclose(fp);
332             sfree(hdr);
333             sfree(suffix);
334             sfree(fnmdup);
335         }
336         /* free pr structure */
337         sfree(prframecurrent->gr);
338         sfree(prframecurrent->r);
339         sfree(prframecurrent);
340         /* free sq structure */
341         sfree(sqframecurrent->q);
342         sfree(sqframecurrent->s);
343         sfree(sqframecurrent);
344     }
345     while (read_next_x(oenv, status, &t, x, box));
346     close_trj(status);
347
348     /* normalize histo */
349     normalize_probability(pr->grn, pr->gr);
350     sq = convert_histogram_to_intensity_curve(pr, start_q, end_q, q_step);
351     /* prepare pr.xvg */
352     fp = xvgropen(opt2fn_null("-pr", NFILE, fnm), "G(r)", "Distance (nm)", "Probability", oenv);
353     for (i = 0; i < pr->grn; i++)
354     {
355         fprintf(fp, "%10.6f%10.6f\n", pr->r[i], pr->gr[i]);
356     }
357     xvgrclose(fp);
358
359     /* prepare sq.xvg */
360     fp = xvgropen(opt2fn_null("-sq", NFILE, fnm), "I(q)", "q (nm^-1)", "s(q)/s(0)", oenv);
361     for (i = 0; i < sq->qn; i++)
362     {
363         fprintf(fp, "%10.6f%10.6f\n", sq->q[i], sq->s[i]);
364     }
365     xvgrclose(fp);
366     /*
367      * Clean up memory
368      */
369     sfree(pr->gr);
370     sfree(pr->r);
371     sfree(pr);
372     sfree(sq->q);
373     sfree(sq->s);
374     sfree(sq);
375
376     please_cite(stdout, "Garmay2012");
377
378     return 0;
379 }