Make PBC type enumeration into PbcType enum class
[gromacs.git] / src / gromacs / gmxana / gmx_filter.cpp
blobc9567a806b3c9a39195219fdc9ebe8ee9b432c2a
1 /*
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,2015,2017,2019,2020, 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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstring>
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/fileio/confio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/gmxana/gmx_ana.h"
46 #include "gromacs/gmxana/princ.h"
47 #include "gromacs/math/do_fit.h"
48 #include "gromacs/math/utilities.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include "gromacs/pbcutil/rmpbc.h"
52 #include "gromacs/topology/index.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/arraysize.h"
55 #include "gromacs/utility/smalloc.h"
57 int gmx_filter(int argc, char* argv[])
59 const char* desc[] = {
60 "[THISMODULE] performs frequency filtering on a trajectory.",
61 "The filter shape is cos([GRK]pi[grk] t/A) + 1 from -A to +A, where A is given",
62 "by the option [TT]-nf[tt] times the time step in the input trajectory.",
63 "This filter reduces fluctuations with period A by 85%, with period",
64 "2*A by 50% and with period 3*A by 17% for low-pass filtering.",
65 "Both a low-pass and high-pass filtered trajectory can be written.[PAR]",
67 "Option [TT]-ol[tt] writes a low-pass filtered trajectory.",
68 "A frame is written every [TT]-nf[tt] input frames.",
69 "This ratio of filter length and output interval ensures a good",
70 "suppression of aliasing of high-frequency motion, which is useful for",
71 "making smooth movies. Also averages of properties which are linear",
72 "in the coordinates are preserved, since all input frames are weighted",
73 "equally in the output.",
74 "When all frames are needed, use the [TT]-all[tt] option.[PAR]",
76 "Option [TT]-oh[tt] writes a high-pass filtered trajectory.",
77 "The high-pass filtered coordinates are added to the coordinates",
78 "from the structure file. When using high-pass filtering use [TT]-fit[tt]",
79 "or make sure you use a trajectory that has been fitted on",
80 "the coordinates in the structure file."
83 static int nf = 10;
84 static gmx_bool bNoJump = TRUE, bFit = FALSE, bLowAll = FALSE;
85 t_pargs pa[] = {
86 { "-nf",
87 FALSE,
88 etINT,
89 { &nf },
90 "Sets the filter length as well as the output interval for low-pass filtering" },
91 { "-all", FALSE, etBOOL, { &bLowAll }, "Write all low-pass filtered frames" },
92 { "-nojump", FALSE, etBOOL, { &bNoJump }, "Remove jumps of atoms across the box" },
93 { "-fit", FALSE, etBOOL, { &bFit }, "Fit all frames to a reference structure" }
95 const char * topfile, *lowfile, *highfile;
96 gmx_bool bTop = FALSE;
97 t_topology top;
98 PbcType pbcType = PbcType::Unset;
99 rvec* xtop;
100 matrix topbox, *box, boxf;
101 char* grpname;
102 int isize;
103 int* index;
104 real* w_rls = nullptr;
105 t_trxstatus* in;
106 t_trxstatus * outl, *outh;
107 int nffr, i, fr, nat, j, d, m;
108 int* ind;
109 real flen, *filt, sum, *t;
110 rvec xcmtop, xcm, **x, *ptr, *xf, *xn, *xp, hbox;
111 gmx_output_env_t* oenv;
112 gmx_rmpbc_t gpbc = nullptr;
114 #define NLEG asize(leg)
115 t_filenm fnm[] = { { efTRX, "-f", nullptr, ffREAD },
116 { efTPS, nullptr, nullptr, ffOPTRD },
117 { efNDX, nullptr, nullptr, ffOPTRD },
118 { efTRO, "-ol", "lowpass", ffOPTWR },
119 { efTRO, "-oh", "highpass", ffOPTWR } };
120 #define NFILE asize(fnm)
122 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW, NFILE, fnm, asize(pa), pa,
123 asize(desc), desc, 0, nullptr, &oenv))
125 return 0;
128 highfile = opt2fn_null("-oh", NFILE, fnm);
129 if (highfile)
131 topfile = ftp2fn(efTPS, NFILE, fnm);
132 lowfile = opt2fn_null("-ol", NFILE, fnm);
134 else
136 topfile = ftp2fn_null(efTPS, NFILE, fnm);
137 lowfile = opt2fn("-ol", NFILE, fnm);
139 if (topfile)
141 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &pbcType, &xtop, nullptr, topbox, TRUE);
142 if (bTop)
144 gpbc = gmx_rmpbc_init(&top.idef, pbcType, top.atoms.nr);
145 gmx_rmpbc(gpbc, top.atoms.nr, topbox, xtop);
149 clear_rvec(xcmtop);
150 if (bFit)
152 fprintf(stderr, "Select group for least squares fit\n");
153 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
154 /* Set the weight */
155 snew(w_rls, top.atoms.nr);
156 for (i = 0; i < isize; i++)
158 w_rls[index[i]] = top.atoms.atom[index[i]].m;
160 calc_xcm(xtop, isize, index, top.atoms.atom, xcmtop, FALSE);
161 for (j = 0; j < top.atoms.nr; j++)
163 rvec_dec(xtop[j], xcmtop);
167 /* The actual filter length flen can actually be any real number */
168 flen = 2 * nf;
169 /* nffr is the number of frames that we filter over */
170 nffr = 2 * nf - 1;
171 snew(filt, nffr);
172 sum = 0;
173 for (i = 0; i < nffr; i++)
175 filt[i] = std::cos(2 * M_PI * (i - nf + 1) / static_cast<real>(flen)) + 1;
176 sum += filt[i];
178 fprintf(stdout, "filter weights:");
179 for (i = 0; i < nffr; i++)
181 filt[i] /= sum;
182 fprintf(stdout, " %5.3f", filt[i]);
184 fprintf(stdout, "\n");
186 snew(t, nffr);
187 snew(x, nffr);
188 snew(box, nffr);
190 nat = read_first_x(oenv, &in, opt2fn("-f", NFILE, fnm), &(t[nffr - 1]), &(x[nffr - 1]), box[nffr - 1]);
191 snew(ind, nat);
192 for (i = 0; i < nat; i++)
194 ind[i] = i;
196 /* x[nffr - 1] was already allocated by read_first_x */
197 for (i = 0; i < nffr - 1; i++)
199 snew(x[i], nat);
201 snew(xf, nat);
202 if (lowfile)
204 outl = open_trx(lowfile, "w");
206 else
208 outl = nullptr;
210 if (highfile)
212 outh = open_trx(highfile, "w");
214 else
216 outh = nullptr;
219 fr = 0;
222 xn = x[nffr - 1];
223 if (bNoJump && fr > 0)
225 xp = x[nffr - 2];
226 for (j = 0; j < nat; j++)
228 for (d = 0; d < DIM; d++)
230 hbox[d] = 0.5 * box[nffr - 1][d][d];
233 for (i = 0; i < nat; i++)
235 for (m = DIM - 1; m >= 0; m--)
237 if (hbox[m] > 0)
239 while (xn[i][m] - xp[i][m] <= -hbox[m])
241 for (d = 0; d <= m; d++)
243 xn[i][d] += box[nffr - 1][m][d];
246 while (xn[i][m] - xp[i][m] > hbox[m])
248 for (d = 0; d <= m; d++)
250 xn[i][d] -= box[nffr - 1][m][d];
257 if (bTop)
259 gmx_rmpbc(gpbc, nat, box[nffr - 1], xn);
261 if (bFit)
263 calc_xcm(xn, isize, index, top.atoms.atom, xcm, FALSE);
264 for (j = 0; j < nat; j++)
266 rvec_dec(xn[j], xcm);
268 do_fit(nat, w_rls, xtop, xn);
269 for (j = 0; j < nat; j++)
271 rvec_inc(xn[j], xcmtop);
274 if (fr >= nffr && (outh || bLowAll || fr % nf == nf - 1))
276 /* Lowpass filtering */
277 for (j = 0; j < nat; j++)
279 clear_rvec(xf[j]);
281 clear_mat(boxf);
282 for (i = 0; i < nffr; i++)
284 for (j = 0; j < nat; j++)
286 for (d = 0; d < DIM; d++)
288 xf[j][d] += filt[i] * x[i][j][d];
291 for (j = 0; j < DIM; j++)
293 for (d = 0; d < DIM; d++)
295 boxf[j][d] += filt[i] * box[i][j][d];
299 if (outl && (bLowAll || fr % nf == nf - 1))
301 write_trx(outl, nat, ind, topfile ? &(top.atoms) : nullptr, 0, t[nf - 1],
302 bFit ? topbox : boxf, xf, nullptr, nullptr);
304 if (outh)
306 /* Highpass filtering */
307 for (j = 0; j < nat; j++)
309 for (d = 0; d < DIM; d++)
311 xf[j][d] = xtop[j][d] + x[nf - 1][j][d] - xf[j][d];
314 if (bFit)
316 for (j = 0; j < nat; j++)
318 rvec_inc(xf[j], xcmtop);
321 for (j = 0; j < DIM; j++)
323 for (d = 0; d < DIM; d++)
325 boxf[j][d] = topbox[j][d] + box[nf - 1][j][d] - boxf[j][d];
328 write_trx(outh, nat, ind, topfile ? &(top.atoms) : nullptr, 0, t[nf - 1],
329 bFit ? topbox : boxf, xf, nullptr, nullptr);
332 /* Cycle all the pointer and the box by one */
333 ptr = x[0];
334 for (i = 0; i < nffr - 1; i++)
336 t[i] = t[i + 1];
337 x[i] = x[i + 1];
338 copy_mat(box[i + 1], box[i]);
340 x[nffr - 1] = ptr;
341 fr++;
342 } while (read_next_x(oenv, in, &(t[nffr - 1]), x[nffr - 1], box[nffr - 1]));
344 if (bTop)
346 gmx_rmpbc_done(gpbc);
349 if (outh)
351 close_trx(outh);
353 if (outl)
355 close_trx(outl);
357 close_trx(in);
359 return 0;