Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_rms.c
index eb83a075df7fc7ff5063b02587ff50f4dae4a67f..dccad9e472458a25f2b947ec369f7cb02bb08ac0 100644 (file)
@@ -1,63 +1,63 @@
 /*
+ * This file is part of the GROMACS molecular simulation package.
  *
- *                This source code is part of
- *
- *                 G   R   O   M   A   C   S
- *
- *          GROningen MAchine for Chemical Simulations
- *
- *                        VERSION 3.2.0
- * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
- * Copyright (c) 2001-2004, The GROMACS development team,
- * check out http://www.gromacs.org for more information.
-
- * This program is free software; you can redistribute it and/or
- * modify it under the terms of the GNU General Public License
- * as published by the Free Software Foundation; either version 2
+ * Copyright (c) 2001-2004, The GROMACS development team.
+ * Copyright (c) 2013,2014, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
  * of the License, or (at your option) any later version.
  *
- * If you want to redistribute modifications, please consider that
- * scientific software is very special. Version control is crucial -
- * bugs must be traceable. We will be happy to consider code for
- * inclusion in the official distribution, but derived work must not
- * be called official GROMACS. Details are found in the README & COPYING
- * files - if they are missing, get the official version at www.gromacs.org.
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
  *
- * To help us fund GROMACS development, we humbly ask that you cite
- * the papers on the package - you can find them in the top README file.
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
  *
- * For more info, check our website at http://www.gromacs.org
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
  *
- * And Hey:
- * Green Red Orange Magenta Azure Cyan Skyblue
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
  */
-#ifdef HAVE_CONFIG_H
-#include <config.h>
-#endif
+#include "gmxpre.h"
 
-#include "smalloc.h"
 #include <math.h>
-#include "macros.h"
-#include "typedefs.h"
-#include "xvgr.h"
-#include "copyrite.h"
-#include "statutil.h"
-#include "string2.h"
-#include "vec.h"
-#include "index.h"
-#include "gmx_fatal.h"
-#include "gromacs/fileio/futil.h"
-#include "princ.h"
-#include "rmpbc.h"
-#include "do_fit.h"
+#include <stdlib.h>
+
+#include "gromacs/commandline/pargs.h"
 #include "gromacs/fileio/matio.h"
 #include "gromacs/fileio/tpxio.h"
 #include "gromacs/fileio/trxio.h"
-#include "cmat.h"
-#include "viewit.h"
-#include "gmx_ana.h"
-
+#include "gromacs/fileio/xvgr.h"
+#include "gromacs/gmxana/cmat.h"
+#include "gromacs/gmxana/gmx_ana.h"
+#include "gromacs/gmxana/princ.h"
+#include "gromacs/legacyheaders/copyrite.h"
+#include "gromacs/legacyheaders/macros.h"
+#include "gromacs/legacyheaders/typedefs.h"
+#include "gromacs/legacyheaders/viewit.h"
+#include "gromacs/math/do_fit.h"
+#include "gromacs/math/vec.h"
+#include "gromacs/pbcutil/rmpbc.h"
+#include "gromacs/topology/index.h"
+#include "gromacs/utility/fatalerror.h"
+#include "gromacs/utility/futil.h"
+#include "gromacs/utility/smalloc.h"
 
 static void norm_princ(t_atoms *atoms, int isize, atom_id *index, int natoms,
                        rvec *x)
@@ -94,10 +94,9 @@ static void norm_princ(t_atoms *atoms, int isize, atom_id *index, int natoms,
 
 int gmx_rms(int argc, char *argv[])
 {
-    const char
-                   *desc[] =
+    const char     *desc[] =
     {
-        "[TT]g_rms[tt] compares two structures by computing the root mean square",
+        "[THISMODULE] compares two structures by computing the root mean square",
         "deviation (RMSD), the size-independent [GRK]rho[grk] similarity parameter",
         "([TT]rho[tt]) or the scaled [GRK]rho[grk] ([TT]rhosc[tt]), ",
         "see Maiorov & Crippen, Proteins [BB]22[bb], 273 (1995).",
@@ -118,7 +117,7 @@ int gmx_rms(int argc, char *argv[])
         "Option [TT]-m[tt] produces a matrix in [TT].xpm[tt] format of",
         "comparison values of each structure in the trajectory with respect to",
         "each other structure. This file can be visualized with for instance",
-        "[TT]xv[tt] and can be converted to postscript with [TT]xpm2ps[tt].[PAR]",
+        "[TT]xv[tt] and can be converted to postscript with [gmx-xpm2ps].[PAR]",
 
         "Option [TT]-fit[tt] controls the least-squares fitting of",
         "the structures on top of each other: complete fit (rotation and",
@@ -228,7 +227,7 @@ int gmx_rms(int argc, char *argv[])
     int             ePBC;
     t_iatom        *iatom = NULL;
 
-    matrix          box;
+    matrix          box = {{0}};
     rvec           *x, *xp, *xm = NULL, **mat_x = NULL, **mat_x2, *mat_x2_j = NULL, vec1,
                     vec2;
     t_trxstatus    *status;
@@ -268,8 +267,8 @@ int gmx_rms(int argc, char *argv[])
     };
 #define NFILE asize(fnm)
 
-    if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW
-                           | PCA_BE_NICE, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL,
+    if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW,
+                           NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL,
                            &oenv))
     {
         return 0;
@@ -322,6 +321,8 @@ int gmx_rms(int argc, char *argv[])
     bPrev = (prev > 0);
     if (bPrev)
     {
+        fprintf(stderr, "WARNING: using option -prev with large trajectories will\n"
+                "         require a lot of memory and could lead to crashes\n");
         prev = abs(prev);
         if (freq != 1)
         {
@@ -1063,11 +1064,11 @@ int gmx_rms(int argc, char *argv[])
                     del_yaxis[i] = delta_maxy*i/del_lev;
                 }
                 sprintf(buf, "%s %s vs. delta t", gn_rms[0], whatname[ewhat]);
-                fp = ffopen("delta.xpm", "w");
+                fp = gmx_ffopen("delta.xpm", "w");
                 write_xpm(fp, 0, buf, "density", output_env_get_time_label(oenv), whatlabel[ewhat],
                           delta_xsize, del_lev+1, del_xaxis, del_yaxis,
                           delta, 0.0, delta_max, rlo, rhi, &nlevels);
-                ffclose(fp);
+                gmx_ffclose(fp);
             }
             if (opt2bSet("-bin", NFILE, fnm))
             {
@@ -1080,7 +1081,7 @@ int gmx_rms(int argc, char *argv[])
                         gmx_fatal(FARGS, "Error writing to output file");
                     }
                 }
-                ffclose(fp);
+                gmx_ffclose(fp);
             }
         }
         if (bBond)
@@ -1135,9 +1136,9 @@ int gmx_rms(int argc, char *argv[])
     for (i = 0; (i < teller); i++)
     {
         if (bSplit && i > 0 &&
-            abs(time[bPrev ? freq*i : i]/output_env_get_time_factor(oenv)) < 1e-5)
+            fabs(time[bPrev ? freq*i : i]/output_env_get_time_factor(oenv)) < 1e-5)
         {
-            fprintf(fp, "&\n");
+            fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
         }
         fprintf(fp, "%12.7f", time[bPrev ? freq*i : i]);
         for (j = 0; (j < nrms); j++)
@@ -1150,7 +1151,7 @@ int gmx_rms(int argc, char *argv[])
         }
         fprintf(fp, "\n");
     }
-    ffclose(fp);
+    gmx_ffclose(fp);
 
     if (bMirror)
     {
@@ -1177,9 +1178,9 @@ int gmx_rms(int argc, char *argv[])
         }
         for (i = 0; (i < teller); i++)
         {
-            if (bSplit && i > 0 && abs(time[i]) < 1e-5)
+            if (bSplit && i > 0 && fabs(time[i]) < 1e-5)
             {
-                fprintf(fp, "&\n");
+                fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             }
             fprintf(fp, "%12.7f", time[i]);
             for (j = 0; (j < nrms); j++)
@@ -1188,7 +1189,7 @@ int gmx_rms(int argc, char *argv[])
             }
             fprintf(fp, "\n");
         }
-        ffclose(fp);
+        gmx_ffclose(fp);
     }
 
     if (bAv)
@@ -1200,7 +1201,7 @@ int gmx_rms(int argc, char *argv[])
         {
             fprintf(fp, "%10d  %10g\n", j, rlstot/teller);
         }
-        ffclose(fp);
+        gmx_ffclose(fp);
     }
 
     if (bNorm)
@@ -1210,7 +1211,7 @@ int gmx_rms(int argc, char *argv[])
         {
             fprintf(fp, "%10d  %10g\n", j, rlsnorm[j]/teller);
         }
-        ffclose(fp);
+        gmx_ffclose(fp);
     }
     do_view(oenv, opt2fn_null("-a", NFILE, fnm), "-graphtype bar");
     do_view(oenv, opt2fn("-o", NFILE, fnm), NULL);