Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_mindist.c
index 9070b9864cae107b79303dfd55b30f28b8e66b3e..920f6cc3a7465e6044e37091d1786049027016b7 100644 (file)
@@ -1,74 +1,90 @@
 /*
+ * 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 <math.h>
 #include <stdlib.h>
-
-#include "sysstuff.h"
 #include <string.h>
-#include "typedefs.h"
-#include "smalloc.h"
-#include "macros.h"
-#include "vec.h"
-#include "xvgr.h"
-#include "pbc.h"
-#include "copyrite.h"
-#include "futil.h"
-#include "statutil.h"
-#include "index.h"
-#include "tpxio.h"
-#include "rmpbc.h"
-#include "xtcio.h"
-#include "gmx_ana.h"
-
-
-static void periodic_dist(matrix box, rvec x[], int n, atom_id index[],
+
+#include "gromacs/commandline/pargs.h"
+#include "gromacs/fileio/tpxio.h"
+#include "gromacs/fileio/trxio.h"
+#include "gromacs/fileio/xvgr.h"
+#include "gromacs/gmxana/gmx_ana.h"
+#include "gromacs/legacyheaders/macros.h"
+#include "gromacs/legacyheaders/names.h"
+#include "gromacs/legacyheaders/typedefs.h"
+#include "gromacs/legacyheaders/viewit.h"
+#include "gromacs/math/vec.h"
+#include "gromacs/pbcutil/pbc.h"
+#include "gromacs/pbcutil/rmpbc.h"
+#include "gromacs/topology/index.h"
+#include "gromacs/utility/futil.h"
+#include "gromacs/utility/smalloc.h"
+
+
+static void periodic_dist(int ePBC,
+                          matrix box, rvec x[], int n, atom_id index[],
                           real *rmin, real *rmax, int *min_ind)
 {
-#define NSHIFT 26
-    int  sx, sy, sz, i, j, s;
+#define NSHIFT_MAX 26
+    int  nsz, nshift, sx, sy, sz, i, j, s;
     real sqr_box, r2min, r2max, r2;
-    rvec shift[NSHIFT], d0, d;
+    rvec shift[NSHIFT_MAX], d0, d;
 
-    sqr_box = sqr(min(box[XX][XX], min(box[YY][YY], box[ZZ][ZZ])));
+    sqr_box = min(norm2(box[XX]), norm2(box[YY]));
+    if (ePBC == epbcXYZ)
+    {
+        sqr_box = min(sqr_box, norm2(box[ZZ]));
+        nsz     = 1;
+    }
+    else if (ePBC == epbcXY)
+    {
+        nsz = 0;
+    }
+    else
+    {
+        gmx_fatal(FARGS, "pbc = %s is not supported by g_mindist",
+                  epbc_names[ePBC]);
+        nsz = 0; /* Keep compilers quiet */
+    }
 
-    s = 0;
-    for (sz = -1; sz <= 1; sz++)
+    nshift = 0;
+    for (sz = -nsz; sz <= nsz; sz++)
     {
         for (sy = -1; sy <= 1; sy++)
         {
@@ -78,9 +94,10 @@ static void periodic_dist(matrix box, rvec x[], int n, atom_id index[],
                 {
                     for (i = 0; i < DIM; i++)
                     {
-                        shift[s][i] = sx*box[XX][i]+sy*box[YY][i]+sz*box[ZZ][i];
+                        shift[nshift][i] =
+                            sx*box[XX][i] + sy*box[YY][i] + sz*box[ZZ][i];
                     }
-                    s++;
+                    nshift++;
                 }
             }
         }
@@ -99,7 +116,7 @@ static void periodic_dist(matrix box, rvec x[], int n, atom_id index[],
             {
                 r2max = r2;
             }
-            for (s = 0; s < NSHIFT; s++)
+            for (s = 0; s < nshift; s++)
             {
                 rvec_add(d0, shift[s], d);
                 r2 = norm2(d);
@@ -150,7 +167,7 @@ static void periodic_mindist_plot(const char *trxfn, const char *outfn,
 
     if (NULL != top)
     {
-        gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms, box);
+        gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
     }
 
     bFirst = TRUE;
@@ -161,7 +178,7 @@ static void periodic_mindist_plot(const char *trxfn, const char *outfn,
             gmx_rmpbc(gpbc, natoms, box, x);
         }
 
-        periodic_dist(box, x, n, index, &rmin, &rmax, ind_min);
+        periodic_dist(ePBC, box, x, n, index, &rmin, &rmax, ind_min);
         if (rmin < rmint)
         {
             rmint    = rmin;
@@ -169,22 +186,22 @@ static void periodic_mindist_plot(const char *trxfn, const char *outfn,
             ind_mini = ind_min[0];
             ind_minj = ind_min[1];
         }
-        if (bSplit && !bFirst && abs(t/output_env_get_time_factor(oenv)) < 1e-5)
+        if (bSplit && !bFirst && fabs(t/output_env_get_time_factor(oenv)) < 1e-5)
         {
-            fprintf(out, "&\n");
+            fprintf(out, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
         }
         fprintf(out, "\t%g\t%6.3f %6.3f %6.3f %6.3f %6.3f\n",
                 output_env_conv_time(oenv, t), rmin, rmax, norm(box[0]), norm(box[1]), norm(box[2]));
         bFirst = FALSE;
     }
-    while (read_next_x(oenv, status, &t, natoms, x, box));
+    while (read_next_x(oenv, status, &t, x, box));
 
     if (NULL != top)
     {
         gmx_rmpbc_done(gpbc);
     }
 
-    ffclose(out);
+    gmx_ffclose(out);
 
     fprintf(stdout,
             "\nThe shortest periodic distance is %g (nm) at time %g (%s),\n"
@@ -335,7 +352,7 @@ void dist_plot(const char *fn, const char *afile, const char *dfile,
     dist = xvgropen(dfile, buf, output_env_get_time_label(oenv), "Distance (nm)", oenv);
     sprintf(buf, "Number of Contacts %s %g nm", bMin ? "<" : ">", rcut);
     num    = nfile ? xvgropen(nfile, buf, output_env_get_time_label(oenv), "Number", oenv) : NULL;
-    atm    = afile ? ffopen(afile, "w") : NULL;
+    atm    = afile ? gmx_ffopen(afile, "w") : NULL;
     trxout = xfile ? open_trx(xfile, "w") : NULL;
 
     if (bMat)
@@ -344,7 +361,7 @@ void dist_plot(const char *fn, const char *afile, const char *dfile,
         {
             snew(leg, 1);
             sprintf(buf, "Internal in %s", grpn[0]);
-            leg[0] = strdup(buf);
+            leg[0] = gmx_strdup(buf);
             xvgr_legend(dist, 0, (const char**)leg, oenv);
             if (num)
             {
@@ -359,7 +376,7 @@ void dist_plot(const char *fn, const char *afile, const char *dfile,
                 for (k = i+1; (k < ng); k++, j++)
                 {
                     sprintf(buf, "%s-%s", grpn[i], grpn[k]);
-                    leg[j] = strdup(buf);
+                    leg[j] = gmx_strdup(buf);
                 }
             }
             xvgr_legend(dist, j, (const char**)leg, oenv);
@@ -375,7 +392,7 @@ void dist_plot(const char *fn, const char *afile, const char *dfile,
         for (i = 0; (i < ng-1); i++)
         {
             sprintf(buf, "%s-%s", grpn[0], grpn[i+1]);
-            leg[i] = strdup(buf);
+            leg[i] = gmx_strdup(buf);
         }
         xvgr_legend(dist, ng-1, (const char**)leg, oenv);
         if (num)
@@ -420,16 +437,16 @@ void dist_plot(const char *fn, const char *afile, const char *dfile,
     bFirst = TRUE;
     do
     {
-        if (bSplit && !bFirst && abs(t/output_env_get_time_factor(oenv)) < 1e-5)
+        if (bSplit && !bFirst && fabs(t/output_env_get_time_factor(oenv)) < 1e-5)
         {
-            fprintf(dist, "&\n");
+            fprintf(dist, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             if (num)
             {
-                fprintf(num, "&\n");
+                fprintf(num, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             }
             if (atm)
             {
-                fprintf(atm, "&\n");
+                fprintf(atm, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
             }
         }
         fprintf(dist, "%12e", output_env_conv_time(oenv, t));
@@ -530,17 +547,17 @@ void dist_plot(const char *fn, const char *afile, const char *dfile,
             fprintf(respertime, "\n");
         }
     }
-    while (read_next_x(oenv, status, &t, natoms, x0, box));
+    while (read_next_x(oenv, status, &t, x0, box));
 
     close_trj(status);
-    ffclose(dist);
+    gmx_ffclose(dist);
     if (num)
     {
-        ffclose(num);
+        gmx_ffclose(num);
     }
     if (atm)
     {
-        ffclose(atm);
+        gmx_ffclose(atm);
     }
     if (trxout)
     {
@@ -599,7 +616,7 @@ int find_residues(t_atoms *atoms, int n, atom_id index[], atom_id **resindex)
     return nres;
 }
 
-void dump_res(FILE *out, int nres, atom_id *resindex, int n, atom_id index[])
+void dump_res(FILE *out, int nres, atom_id *resindex, atom_id index[])
 {
     int i, j;
 
@@ -617,7 +634,7 @@ void dump_res(FILE *out, int nres, atom_id *resindex, int n, atom_id index[])
 int gmx_mindist(int argc, char *argv[])
 {
     const char     *desc[] = {
-        "[TT]g_mindist[tt] computes the distance between one group and a number of",
+        "[THISMODULE] computes the distance between one group and a number of",
         "other groups. Both the minimum distance",
         "(between any pair of atoms from the respective groups)",
         "and the number of contacts within a given",
@@ -633,8 +650,7 @@ int gmx_mindist(int argc, char *argv[])
         "each direction is considered, giving a total of 26 shifts.",
         "It also plots the maximum distance within the group and the lengths",
         "of the three box vectors.[PAR]",
-        "Other programs that calculate distances are [TT]g_dist[tt]",
-        "and [TT]g_bond[tt]."
+        "Also [gmx-distance] calculates distances."
     };
     const char     *bugs[] = {
         "The [TT]-pi[tt] option is very slow."
@@ -694,9 +710,12 @@ int gmx_mindist(int argc, char *argv[])
     };
 #define NFILE asize(fnm)
 
-    parse_common_args(&argc, argv,
-                      PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE,
-                      NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
+    if (!parse_common_args(&argc, argv,
+                           PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT,
+                           NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
+    {
+        return 0;
+    }
 
     trxfnm  = ftp2fn(efTRX, NFILE, fnm);
     ndxfnm  = ftp2fn_null(efNDX, NFILE, fnm);
@@ -769,7 +788,7 @@ int gmx_mindist(int argc, char *argv[])
                              gnx[0], index[0], &residues);
         if (debug)
         {
-            dump_res(debug, nres, residues, gnx[0], index[0]);
+            dump_res(debug, nres, residues, index[0]);
         }
     }
 
@@ -791,7 +810,5 @@ int gmx_mindist(int argc, char *argv[])
         do_view(oenv, numfnm, "-nxy");
     }
 
-    thanx(stderr);
-
     return 0;
 }