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