2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/legacyheaders/typedefs.h"
43 #include "gromacs/utility/smalloc.h"
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/topology/index.h"
47 #include "gromacs/fileio/tpxio.h"
48 #include "gromacs/fileio/trxio.h"
50 #include "gromacs/pbcutil/rmpbc.h"
53 #include "gromacs/math/do_fit.h"
55 int gmx_filter(int argc, char *argv[])
57 const char *desc[] = {
58 "[THISMODULE] performs frequency filtering on a trajectory.",
59 "The filter shape is cos([GRK]pi[grk] t/A) + 1 from -A to +A, where A is given",
60 "by the option [TT]-nf[tt] times the time step in the input trajectory.",
61 "This filter reduces fluctuations with period A by 85%, with period",
62 "2*A by 50% and with period 3*A by 17% for low-pass filtering.",
63 "Both a low-pass and high-pass filtered trajectory can be written.[PAR]",
65 "Option [TT]-ol[tt] writes a low-pass filtered trajectory.",
66 "A frame is written every [TT]-nf[tt] input frames.",
67 "This ratio of filter length and output interval ensures a good",
68 "suppression of aliasing of high-frequency motion, which is useful for",
69 "making smooth movies. Also averages of properties which are linear",
70 "in the coordinates are preserved, since all input frames are weighted",
71 "equally in the output.",
72 "When all frames are needed, use the [TT]-all[tt] option.[PAR]",
74 "Option [TT]-oh[tt] writes a high-pass filtered trajectory.",
75 "The high-pass filtered coordinates are added to the coordinates",
76 "from the structure file. When using high-pass filtering use [TT]-fit[tt]",
77 "or make sure you use a trajectory that has been fitted on",
78 "the coordinates in the structure file."
82 static gmx_bool bNoJump = TRUE, bFit = FALSE, bLowAll = FALSE;
84 { "-nf", FALSE, etINT, {&nf},
85 "Sets the filter length as well as the output interval for low-pass filtering" },
86 { "-all", FALSE, etBOOL, {&bLowAll},
87 "Write all low-pass filtered frames" },
88 { "-nojump", FALSE, etBOOL, {&bNoJump},
89 "Remove jumps of atoms across the box" },
90 { "-fit", FALSE, etBOOL, {&bFit},
91 "Fit all frames to a reference structure" }
93 const char *topfile, *lowfile, *highfile;
94 gmx_bool bTop = FALSE;
98 matrix topbox, *box, boxf;
99 char title[256], *grpname;
104 t_trxstatus *outl, *outh;
105 int nffr, i, fr, nat, j, d, m;
107 real flen, *filt, sum, *t;
108 rvec xcmtop, xcm, **x, *ptr, *xf, *xn, *xp, hbox;
110 gmx_rmpbc_t gpbc = NULL;
112 #define NLEG asize(leg)
114 { efTRX, "-f", NULL, ffREAD },
115 { efTPS, NULL, NULL, ffOPTRD },
116 { efNDX, NULL, NULL, ffOPTRD },
117 { efTRO, "-ol", "lowpass", ffOPTWR },
118 { efTRO, "-oh", "highpass", ffOPTWR }
120 #define NFILE asize(fnm)
122 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
123 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
128 highfile = opt2fn_null("-oh", NFILE, fnm);
131 topfile = ftp2fn(efTPS, NFILE, fnm);
132 lowfile = opt2fn_null("-ol", NFILE, fnm);
136 topfile = ftp2fn_null(efTPS, NFILE, fnm);
137 lowfile = opt2fn("-ol", NFILE, fnm);
141 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC,
142 &xtop, NULL, topbox, TRUE);
145 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
146 gmx_rmpbc(gpbc, top.atoms.nr, topbox, xtop);
153 fprintf(stderr, "Select group for least squares fit\n");
154 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
156 snew(w_rls, top.atoms.nr);
157 for (i = 0; i < isize; i++)
159 w_rls[index[i]] = top.atoms.atom[index[i]].m;
161 calc_xcm(xtop, isize, index, top.atoms.atom, xcmtop, FALSE);
162 for (j = 0; j < top.atoms.nr; j++)
164 rvec_dec(xtop[j], xcmtop);
168 /* The actual filter length flen can actually be any real number */
170 /* nffr is the number of frames that we filter over */
174 for (i = 0; i < nffr; i++)
176 filt[i] = cos(2*M_PI*(i - nf + 1)/(real)flen) + 1;
179 fprintf(stdout, "filter weights:");
180 for (i = 0; i < nffr; i++)
183 fprintf(stdout, " %5.3f", filt[i]);
185 fprintf(stdout, "\n");
191 nat = read_first_x(oenv, &in, opt2fn("-f", NFILE, fnm),
192 &(t[nffr - 1]), &(x[nffr - 1]), box[nffr - 1]);
194 for (i = 0; i < nat; i++)
198 /* x[nffr - 1] was already allocated by read_first_x */
199 for (i = 0; i < nffr-1; i++)
206 outl = open_trx(lowfile, "w");
214 outh = open_trx(highfile, "w");
225 if (bNoJump && fr > 0)
228 for (j = 0; j < nat; j++)
230 for (d = 0; d < DIM; d++)
232 hbox[d] = 0.5*box[nffr - 1][d][d];
235 for (i = 0; i < nat; i++)
237 for (m = DIM-1; m >= 0; m--)
241 while (xn[i][m] - xp[i][m] <= -hbox[m])
243 for (d = 0; d <= m; d++)
245 xn[i][d] += box[nffr - 1][m][d];
248 while (xn[i][m] - xp[i][m] > hbox[m])
250 for (d = 0; d <= m; d++)
252 xn[i][d] -= box[nffr - 1][m][d];
261 gmx_rmpbc(gpbc, nat, box[nffr - 1], xn);
265 calc_xcm(xn, isize, index, top.atoms.atom, xcm, FALSE);
266 for (j = 0; j < nat; j++)
268 rvec_dec(xn[j], xcm);
270 do_fit(nat, w_rls, xtop, xn);
271 for (j = 0; j < nat; j++)
273 rvec_inc(xn[j], xcmtop);
276 if (fr >= nffr && (outh || bLowAll || fr % nf == nf - 1))
278 /* Lowpass filtering */
279 for (j = 0; j < nat; j++)
284 for (i = 0; i < nffr; i++)
286 for (j = 0; j < nat; j++)
288 for (d = 0; d < DIM; d++)
290 xf[j][d] += filt[i]*x[i][j][d];
293 for (j = 0; j < DIM; j++)
295 for (d = 0; d < DIM; d++)
297 boxf[j][d] += filt[i]*box[i][j][d];
301 if (outl && (bLowAll || fr % nf == nf - 1))
303 write_trx(outl, nat, ind, topfile ? &(top.atoms) : NULL,
304 0, t[nf - 1], bFit ? topbox : boxf, xf, NULL, NULL);
308 /* Highpass filtering */
309 for (j = 0; j < nat; j++)
311 for (d = 0; d < DIM; d++)
313 xf[j][d] = xtop[j][d] + x[nf - 1][j][d] - xf[j][d];
318 for (j = 0; j < nat; j++)
320 rvec_inc(xf[j], xcmtop);
323 for (j = 0; j < DIM; j++)
325 for (d = 0; d < DIM; d++)
327 boxf[j][d] = topbox[j][d] + box[nf - 1][j][d] - boxf[j][d];
330 write_trx(outh, nat, ind, topfile ? &(top.atoms) : NULL,
331 0, t[nf - 1], bFit ? topbox : boxf, xf, NULL, NULL);
334 /* Cycle all the pointer and the box by one */
336 for (i = 0; i < nffr-1; i++)
340 copy_mat(box[i+1], box[i]);
345 while (read_next_x(oenv, in, &(t[nffr - 1]), x[nffr - 1], box[nffr - 1]));
349 gmx_rmpbc_done(gpbc);