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