2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
43 #include "gromacs/fileio/matio.h"
44 #include "gromacs/fileio/xvgr.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/utility/futil.h"
47 #include "gromacs/utility/smalloc.h"
49 t_mat
*init_mat(int n1
, gmx_bool b1D
)
60 m
->mat
= mk_matrix(n1
, n1
, b1D
);
69 void copy_t_mat(t_mat
*dst
, t_mat
*src
)
73 if (dst
->nn
!= src
->nn
)
75 fprintf(stderr
, "t_mat structures not identical in size dst %d src %d\n", dst
->nn
, src
->nn
);
78 dst
->maxrms
= src
->maxrms
;
79 dst
->minrms
= src
->minrms
;
80 dst
->sumrms
= src
->sumrms
;
81 for (i
= 0; (i
< src
->nn
); i
++)
83 for (j
= 0; (j
< src
->nn
); j
++)
85 dst
->mat
[i
][j
] = src
->mat
[i
][j
];
87 dst
->erow
[i
] = src
->erow
[i
];
88 dst
->m_ind
[i
] = src
->m_ind
[i
];
92 void enlarge_mat(t_mat
*m
, int deltan
)
96 srenew(m
->erow
, m
->nn
+deltan
);
97 srenew(m
->m_ind
, m
->nn
+deltan
);
98 srenew(m
->mat
, m
->nn
+deltan
);
100 /* Reallocate existing rows in the matrix, and set them to zero */
101 for (i
= 0; (i
< m
->nn
); i
++)
103 srenew(m
->mat
[i
], m
->nn
+deltan
);
104 for (j
= m
->nn
; (j
< m
->nn
+deltan
); j
++)
109 /* Allocate new rows of the matrix, set energies to zero */
110 for (i
= m
->nn
; (i
< m
->nn
+deltan
); i
++)
113 snew(m
->mat
[i
], m
->nn
+deltan
);
118 void reset_index(t_mat
*m
)
122 for (i
= 0; (i
< m
->n1
); i
++)
128 void set_mat_entry(t_mat
*m
, int i
, int j
, real val
)
130 m
->mat
[i
][j
] = m
->mat
[j
][i
] = val
;
131 m
->maxrms
= std::max(m
->maxrms
, val
);
134 m
->minrms
= std::min(m
->minrms
, val
);
137 m
->nn
= std::max(m
->nn
, std::max(j
+1, i
+1));
140 void done_mat(t_mat
**m
)
142 done_matrix((*m
)->n1
, &((*m
)->mat
));
149 real
mat_energy(t_mat
*m
)
154 for (j
= 0; (j
< m
->nn
-1); j
++)
156 emat
+= sqr(m
->mat
[j
][j
+1]);
161 void swap_rows(t_mat
*m
, int iswap
, int jswap
)
167 itmp
= m
->m_ind
[iswap
];
168 m
->m_ind
[iswap
] = m
->m_ind
[jswap
];
169 m
->m_ind
[jswap
] = itmp
;
171 /* Swap rows (since the matrix is an array of pointers) */
173 m
->mat
[iswap
] = m
->mat
[jswap
];
177 for (i
= 0; (i
< m
->nn
); i
++)
179 ttt
= m
->mat
[i
][iswap
];
180 m
->mat
[i
][iswap
] = m
->mat
[i
][jswap
];
181 m
->mat
[i
][jswap
] = ttt
;
185 void swap_mat(t_mat
*m
)
190 tmp
= init_mat(m
->nn
, FALSE
);
191 for (i
= 0; (i
< m
->nn
); i
++)
193 for (j
= 0; (j
< m
->nn
); j
++)
195 tmp
->mat
[m
->m_ind
[i
]][m
->m_ind
[j
]] = m
->mat
[i
][j
];
198 /*tmp->mat[i][j] = m->mat[m->m_ind[i]][m->m_ind[j]]; */
199 for (i
= 0; (i
< m
->nn
); i
++)
201 for (j
= 0; (j
< m
->nn
); j
++)
203 m
->mat
[i
][j
] = tmp
->mat
[i
][j
];
209 void low_rmsd_dist(const char *fn
, real maxrms
, int nn
, real
**mat
,
210 const gmx_output_env_t
*oenv
)
218 for (i
= 0; i
< nn
; i
++)
220 for (j
= i
+1; j
< nn
; j
++)
222 x
= static_cast<int>(fac
*mat
[i
][j
]+0.5);
230 fp
= xvgropen(fn
, "RMS Distribution", "RMS (nm)", "a.u.", oenv
);
231 for (i
= 0; (i
< 101); i
++)
233 fprintf(fp
, "%10g %10d\n", i
/fac
, histo
[i
]);
239 void rmsd_distribution(const char *fn
, t_mat
*rms
, const gmx_output_env_t
*oenv
)
241 low_rmsd_dist(fn
, rms
->maxrms
, rms
->nn
, rms
->mat
, oenv
);
244 t_clustid
*new_clustid(int n1
)
250 for (i
= 0; (i
< n1
); i
++)