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