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