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