b44e42f979a497ab2425d8596dbaaca39c773eac
[alexxy/gromacs.git] / src / gromacs / math / vecdump.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2016, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "vecdump.h"
40
41 #include <cstdio>
42 #include <cstdlib>
43
44 #include "gromacs/math/vec.h"
45 #include "gromacs/utility/stringutil.h"
46 #include "gromacs/utility/txtdump.h"
47
48 void pr_ivec(FILE *fp, int indent, const char *title, const int vec[], int n, gmx_bool bShowNumbers)
49 {
50     int i;
51
52     if (available(fp, vec, indent, title))
53     {
54         indent = pr_title_n(fp, indent, title, n);
55         for (i = 0; i < n; i++)
56         {
57             pr_indent(fp, indent);
58             fprintf(fp, "%s[%d]=%d\n", title, bShowNumbers ? i : -1, vec[i]);
59         }
60     }
61 }
62
63 void pr_ivec_block(FILE *fp, int indent, const char *title, const int vec[], int n, gmx_bool bShowNumbers)
64 {
65     int i, j;
66
67     if (available(fp, vec, indent, title))
68     {
69         indent = pr_title_n(fp, indent, title, n);
70         i      = 0;
71         while (i < n)
72         {
73             j = i+1;
74             while (j < n && vec[j] == vec[j-1]+1)
75             {
76                 j++;
77             }
78             /* Print consecutive groups of 3 or more as blocks */
79             if (j - i < 3)
80             {
81                 while (i < j)
82                 {
83                     pr_indent(fp, indent);
84                     fprintf(fp, "%s[%d]=%d\n",
85                             title, bShowNumbers ? i : -1, vec[i]);
86                     i++;
87                 }
88             }
89             else
90             {
91                 pr_indent(fp, indent);
92                 fprintf(fp, "%s[%d,...,%d] = {%d,...,%d}\n",
93                         title,
94                         bShowNumbers ? i : -1,
95                         bShowNumbers ? j-1 : -1,
96                         vec[i], vec[j-1]);
97                 i = j;
98             }
99         }
100     }
101 }
102
103 void pr_bvec(FILE *fp, int indent, const char *title, const gmx_bool vec[], int n, gmx_bool bShowNumbers)
104 {
105     int i;
106
107     if (available(fp, vec, indent, title))
108     {
109         indent = pr_title_n(fp, indent, title, n);
110         for (i = 0; i < n; i++)
111         {
112             pr_indent(fp, indent);
113             fprintf(fp, "%s[%d]=%s\n", title, bShowNumbers ? i : -1,
114                     gmx::boolToString(vec[i]));
115         }
116     }
117 }
118
119 void pr_ivecs(FILE *fp, int indent, const char *title, const ivec vec[], int n, gmx_bool bShowNumbers)
120 {
121     int i, j;
122
123     if (available(fp, vec, indent, title))
124     {
125         indent = pr_title_nxn(fp, indent, title, n, DIM);
126         for (i = 0; i < n; i++)
127         {
128             pr_indent(fp, indent);
129             fprintf(fp, "%s[%d]={", title, bShowNumbers ? i : -1);
130             for (j = 0; j < DIM; j++)
131             {
132                 if (j != 0)
133                 {
134                     fprintf(fp, ", ");
135                 }
136                 fprintf(fp, "%d", vec[i][j]);
137             }
138             fprintf(fp, "}\n");
139         }
140     }
141 }
142
143 template <typename T>
144 static void printRealVector(FILE *fp, int indent, const char *title, const T vec[], int n, gmx_bool bShowNumbers)
145 {
146     if (available(fp, vec, indent, title))
147     {
148         indent = pr_title_n(fp, indent, title, n);
149         for (int i = 0; i < n; i++)
150         {
151             pr_indent(fp, indent);
152             fprintf(fp, "%s[%d]=%12.5e\n", title, bShowNumbers ? i : -1, vec[i]);
153         }
154     }
155 }
156
157 void pr_rvec(FILE *fp, int indent, const char *title, const real vec[], int n, gmx_bool bShowNumbers)
158 {
159     printRealVector<real>(fp, indent, title, vec, n, bShowNumbers);
160 }
161
162 void pr_fvec(FILE *fp, int indent, const char *title, const float vec[], int n, gmx_bool bShowNumbers)
163 {
164     printRealVector<float>(fp, indent, title, vec, n, bShowNumbers);
165 }
166
167 void pr_dvec(FILE *fp, int indent, const char *title, const double vec[], int n, gmx_bool bShowNumbers)
168 {
169     printRealVector<double>(fp, indent, title, vec, n, bShowNumbers);
170 }
171
172
173 void pr_rvecs_len(FILE *fp, int indent, const char *title, const rvec vec[], int n)
174 {
175     int i, j;
176
177     if (available(fp, vec, indent, title))
178     {
179         indent = pr_title_nxn(fp, indent, title, n, DIM);
180         for (i = 0; i < n; i++)
181         {
182             pr_indent(fp, indent);
183             fprintf(fp, "%s[%5d]={", title, i);
184             for (j = 0; j < DIM; j++)
185             {
186                 if (j != 0)
187                 {
188                     fprintf(fp, ", ");
189                 }
190                 fprintf(fp, "%12.5e", vec[i][j]);
191             }
192             fprintf(fp, "} len=%12.5e\n", norm(vec[i]));
193         }
194     }
195 }
196
197 void pr_rvecs(FILE *fp, int indent, const char *title, const rvec vec[], int n)
198 {
199     const char *fshort = "%12.5e";
200     const char *flong  = "%15.8e";
201     const char *format;
202     int         i, j;
203
204     if (getenv("GMX_PRINT_LONGFORMAT") != NULL)
205     {
206         format = flong;
207     }
208     else
209     {
210         format = fshort;
211     }
212
213     if (available(fp, vec, indent, title))
214     {
215         indent = pr_title_nxn(fp, indent, title, n, DIM);
216         for (i = 0; i < n; i++)
217         {
218             pr_indent(fp, indent);
219             fprintf(fp, "%s[%5d]={", title, i);
220             for (j = 0; j < DIM; j++)
221             {
222                 if (j != 0)
223                 {
224                     fprintf(fp, ", ");
225                 }
226                 fprintf(fp, format, vec[i][j]);
227             }
228             fprintf(fp, "}\n");
229         }
230     }
231 }
232
233
234 void pr_rvecs_of_dim(FILE *fp, int indent, const char *title, const rvec vec[], int n, int dim)
235 {
236     const char *fshort = "%12.5e";
237     const char *flong  = "%15.8e";
238     const char *format;
239     int         i, j;
240
241     if (getenv("GMX_PRINT_LONGFORMAT") != NULL)
242     {
243         format = flong;
244     }
245     else
246     {
247         format = fshort;
248     }
249
250     if (available(fp, vec, indent, title))
251     {
252         indent = pr_title_nxn(fp, indent, title, n, dim);
253         for (i = 0; i < n; i++)
254         {
255             pr_indent(fp, indent);
256             fprintf(fp, "%s[%5d]={", title, i);
257             for (j = 0; j < dim; j++)
258             {
259                 if (j != 0)
260                 {
261                     fprintf(fp, ", ");
262                 }
263                 fprintf(fp, format, vec[i][j]);
264             }
265             fprintf(fp, "}\n");
266         }
267     }
268 }