Merge branch release-4-6
[alexxy/gromacs.git] / src / gromacs / gmxana / cmat.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "cmat.h"
40 #include "smalloc.h"
41 #include "macros.h"
42 #include "vec.h"
43 #include "xvgr.h"
44 #include "matio.h"
45 #include "futil.h"
46
47 t_mat *init_mat(int n1, gmx_bool b1D)
48 {
49     t_mat *m;
50
51     snew(m, 1);
52     m->n1     = n1;
53     m->nn     = 0;
54     m->b1D    = b1D;
55     m->maxrms = 0;
56     m->minrms = 1e20;
57     m->sumrms = 0;
58     m->mat    = mk_matrix(n1, n1, b1D);
59
60     snew(m->erow, n1);
61     snew(m->m_ind, n1);
62     reset_index(m);
63
64     return m;
65 }
66
67 void copy_t_mat(t_mat *dst, t_mat *src)
68 {
69     int i, j;
70
71     if (dst->nn != src->nn)
72     {
73         fprintf(stderr, "t_mat structures not identical in size dst %d src %d\n",dst->nn,src->nn);
74         return;
75     }
76     dst->maxrms = src->maxrms;
77     dst->minrms = src->minrms;
78     dst->sumrms = src->sumrms;
79     for(i = 0; (i < src->nn); i++)
80     {
81         for(j = 0; (j < src->nn); j++)
82         {
83             dst->mat[i][j] = src->mat[i][j];
84         }
85         dst->erow[i] = src->erow[i];
86         dst->m_ind[i] = src->m_ind[i];
87     }
88 }
89
90 void enlarge_mat(t_mat *m, int deltan)
91 {
92     int i, j;
93
94     srenew(m->erow, m->nn+deltan);
95     srenew(m->m_ind, m->nn+deltan);
96     srenew(m->mat, m->nn+deltan);
97
98     /* Reallocate existing rows in the matrix, and set them to zero */
99     for (i = 0; (i < m->nn); i++)
100     {
101         srenew(m->mat[i], m->nn+deltan);
102         for (j = m->nn; (j < m->nn+deltan); j++)
103         {
104             m->mat[i][j] = 0;
105         }
106     }
107     /* Allocate new rows of the matrix, set energies to zero */
108     for (i = m->nn; (i < m->nn+deltan); i++)
109     {
110         m->erow[i]  = 0;
111         snew(m->mat[i], m->nn+deltan);
112     }
113     m->nn += deltan;
114 }
115
116 void reset_index(t_mat *m)
117 {
118     int i;
119
120     for (i = 0; (i < m->n1); i++)
121     {
122         m->m_ind[i] = i;
123     }
124 }
125
126 void set_mat_entry(t_mat *m, int i, int j, real val)
127 {
128     m->mat[i][j] = m->mat[j][i] = val;
129     m->maxrms    = max(m->maxrms, val);
130     if (j != i)
131     {
132         m->minrms  = min(m->minrms, val);
133     }
134     m->sumrms   += val;
135     m->nn        = max(m->nn, max(j+1, i+1));
136 }
137
138 void done_mat(t_mat **m)
139 {
140     done_matrix((*m)->n1, &((*m)->mat));
141     sfree((*m)->m_ind);
142     sfree((*m)->erow);
143     sfree(*m);
144     *m = NULL;
145 }
146
147 real mat_energy(t_mat *m)
148 {
149     int  j;
150     real emat = 0;
151
152     for (j = 0; (j < m->nn-1); j++)
153     {
154         emat += sqr(m->mat[j][j+1]);
155     }
156     return emat;
157 }
158
159 void swap_rows(t_mat *m, int iswap, int jswap)
160 {
161     real *tmp, ttt;
162     int   i, itmp;
163
164     /* Swap indices */
165     itmp            = m->m_ind[iswap];
166     m->m_ind[iswap] = m->m_ind[jswap];
167     m->m_ind[jswap] = itmp;
168
169     /* Swap rows (since the matrix is an array of pointers) */
170     tmp           = m->mat[iswap];
171     m->mat[iswap] = m->mat[jswap];
172     m->mat[jswap] = tmp;
173
174     /* Swap columns */
175     for (i = 0; (i < m->nn); i++)
176     {
177         ttt              = m->mat[i][iswap];
178         m->mat[i][iswap] = m->mat[i][jswap];
179         m->mat[i][jswap] = ttt;
180     }
181 }
182
183 void swap_mat(t_mat *m)
184 {
185     t_mat *tmp;
186     int    i, j;
187
188     tmp = init_mat(m->nn, FALSE);
189     for (i = 0; (i < m->nn); i++)
190     {
191         for (j = 0; (j < m->nn); j++)
192         {
193             tmp->mat[m->m_ind[i]][m->m_ind[j]] = m->mat[i][j];
194         }
195     }
196     /*tmp->mat[i][j] =  m->mat[m->m_ind[i]][m->m_ind[j]]; */
197     for (i = 0; (i < m->nn); i++)
198     {
199         for (j = 0; (j < m->nn); j++)
200         {
201             m->mat[i][j] = tmp->mat[i][j];
202         }
203     }
204     done_mat(&tmp);
205 }
206
207 void low_rmsd_dist(const char *fn, real maxrms, int nn, real **mat,
208                    const output_env_t oenv)
209 {
210     FILE   *fp;
211     int     i, j, *histo, x;
212     real    fac;
213
214     fac = 100/maxrms;
215     snew(histo, 101);
216     for (i = 0; i < nn; i++)
217     {
218         for (j = i+1; j < nn; j++)
219         {
220             x = (int)(fac*mat[i][j]+0.5);
221             if (x <= 100)
222             {
223                 histo[x]++;
224             }
225         }
226     }
227
228     fp = xvgropen(fn, "RMS Distribution", "RMS (nm)", "a.u.", oenv);
229     for (i = 0; (i < 101); i++)
230     {
231         fprintf(fp, "%10g  %10d\n", i/fac, histo[i]);
232     }
233     ffclose(fp);
234     sfree(histo);
235 }
236
237 void rmsd_distribution(const char *fn, t_mat *rms, const output_env_t oenv)
238 {
239     low_rmsd_dist(fn, rms->maxrms, rms->nn, rms->mat, oenv);
240 }
241
242 t_clustid *new_clustid(int n1)
243 {
244     t_clustid *c;
245     int        i;
246
247     snew(c, n1);
248     for (i = 0; (i < n1); i++)
249     {
250         c[i].conf  = i;
251         c[i].clust = i;
252     }
253     return c;
254 }