Moved timing source to C++
[alexxy/gromacs.git] / src / gromacs / gmxana / cmat.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, 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 "cmat.h"
40
41 #include <algorithm>
42
43 #include "gromacs/fileio/matio.h"
44 #include "gromacs/fileio/xvgr.h"
45 #include "gromacs/legacyheaders/macros.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/utility/futil.h"
48 #include "gromacs/utility/smalloc.h"
49
50 t_mat *init_mat(int n1, gmx_bool b1D)
51 {
52     t_mat *m;
53
54     snew(m, 1);
55     m->n1     = n1;
56     m->nn     = 0;
57     m->b1D    = b1D;
58     m->maxrms = 0;
59     m->minrms = 1e20;
60     m->sumrms = 0;
61     m->mat    = mk_matrix(n1, n1, b1D);
62
63     snew(m->erow, n1);
64     snew(m->m_ind, n1);
65     reset_index(m);
66
67     return m;
68 }
69
70 void copy_t_mat(t_mat *dst, t_mat *src)
71 {
72     int i, j;
73
74     if (dst->nn != src->nn)
75     {
76         fprintf(stderr, "t_mat structures not identical in size dst %d src %d\n", dst->nn, src->nn);
77         return;
78     }
79     dst->maxrms = src->maxrms;
80     dst->minrms = src->minrms;
81     dst->sumrms = src->sumrms;
82     for (i = 0; (i < src->nn); i++)
83     {
84         for (j = 0; (j < src->nn); j++)
85         {
86             dst->mat[i][j] = src->mat[i][j];
87         }
88         dst->erow[i]  = src->erow[i];
89         dst->m_ind[i] = src->m_ind[i];
90     }
91 }
92
93 void enlarge_mat(t_mat *m, int deltan)
94 {
95     int i, j;
96
97     srenew(m->erow, m->nn+deltan);
98     srenew(m->m_ind, m->nn+deltan);
99     srenew(m->mat, m->nn+deltan);
100
101     /* Reallocate existing rows in the matrix, and set them to zero */
102     for (i = 0; (i < m->nn); i++)
103     {
104         srenew(m->mat[i], m->nn+deltan);
105         for (j = m->nn; (j < m->nn+deltan); j++)
106         {
107             m->mat[i][j] = 0;
108         }
109     }
110     /* Allocate new rows of the matrix, set energies to zero */
111     for (i = m->nn; (i < m->nn+deltan); i++)
112     {
113         m->erow[i]  = 0;
114         snew(m->mat[i], m->nn+deltan);
115     }
116     m->nn += deltan;
117 }
118
119 void reset_index(t_mat *m)
120 {
121     int i;
122
123     for (i = 0; (i < m->n1); i++)
124     {
125         m->m_ind[i] = i;
126     }
127 }
128
129 void set_mat_entry(t_mat *m, int i, int j, real val)
130 {
131     m->mat[i][j] = m->mat[j][i] = val;
132     m->maxrms    = std::max(m->maxrms, val);
133     if (j != i)
134     {
135         m->minrms  = std::min(m->minrms, val);
136     }
137     m->sumrms   += val;
138     m->nn        = std::max(m->nn, std::max(j+1, i+1));
139 }
140
141 void done_mat(t_mat **m)
142 {
143     done_matrix((*m)->n1, &((*m)->mat));
144     sfree((*m)->m_ind);
145     sfree((*m)->erow);
146     sfree(*m);
147     *m = NULL;
148 }
149
150 real mat_energy(t_mat *m)
151 {
152     int  j;
153     real emat = 0;
154
155     for (j = 0; (j < m->nn-1); j++)
156     {
157         emat += sqr(m->mat[j][j+1]);
158     }
159     return emat;
160 }
161
162 void swap_rows(t_mat *m, int iswap, int jswap)
163 {
164     real *tmp, ttt;
165     int   i, itmp;
166
167     /* Swap indices */
168     itmp            = m->m_ind[iswap];
169     m->m_ind[iswap] = m->m_ind[jswap];
170     m->m_ind[jswap] = itmp;
171
172     /* Swap rows (since the matrix is an array of pointers) */
173     tmp           = m->mat[iswap];
174     m->mat[iswap] = m->mat[jswap];
175     m->mat[jswap] = tmp;
176
177     /* Swap columns */
178     for (i = 0; (i < m->nn); i++)
179     {
180         ttt              = m->mat[i][iswap];
181         m->mat[i][iswap] = m->mat[i][jswap];
182         m->mat[i][jswap] = ttt;
183     }
184 }
185
186 void swap_mat(t_mat *m)
187 {
188     t_mat *tmp;
189     int    i, j;
190
191     tmp = init_mat(m->nn, FALSE);
192     for (i = 0; (i < m->nn); i++)
193     {
194         for (j = 0; (j < m->nn); j++)
195         {
196             tmp->mat[m->m_ind[i]][m->m_ind[j]] = m->mat[i][j];
197         }
198     }
199     /*tmp->mat[i][j] =  m->mat[m->m_ind[i]][m->m_ind[j]]; */
200     for (i = 0; (i < m->nn); i++)
201     {
202         for (j = 0; (j < m->nn); j++)
203         {
204             m->mat[i][j] = tmp->mat[i][j];
205         }
206     }
207     done_mat(&tmp);
208 }
209
210 void low_rmsd_dist(const char *fn, real maxrms, int nn, real **mat,
211                    const output_env_t oenv)
212 {
213     FILE   *fp;
214     int     i, j, *histo, x;
215     real    fac;
216
217     fac = 100/maxrms;
218     snew(histo, 101);
219     for (i = 0; i < nn; i++)
220     {
221         for (j = i+1; j < nn; j++)
222         {
223             x = (int)(fac*mat[i][j]+0.5);
224             if (x <= 100)
225             {
226                 histo[x]++;
227             }
228         }
229     }
230
231     fp = xvgropen(fn, "RMS Distribution", "RMS (nm)", "a.u.", oenv);
232     for (i = 0; (i < 101); i++)
233     {
234         fprintf(fp, "%10g  %10d\n", i/fac, histo[i]);
235     }
236     xvgrclose(fp);
237     sfree(histo);
238 }
239
240 void rmsd_distribution(const char *fn, t_mat *rms, const output_env_t oenv)
241 {
242     low_rmsd_dist(fn, rms->maxrms, rms->nn, rms->mat, oenv);
243 }
244
245 t_clustid *new_clustid(int n1)
246 {
247     t_clustid *c;
248     int        i;
249
250     snew(c, n1);
251     for (i = 0; (i < n1); i++)
252     {
253         c[i].conf  = i;
254         c[i].clust = i;
255     }
256     return c;
257 }