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) 2012,2014,2015,2017,2018 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 /* This file is completely threadsafe - keep it that way! */
45 #include "gromacs/utility/fatalerror.h"
46 #include "gromacs/utility/smalloc.h"
48 static inline void do_rotate(double** a
, int i
, int j
, int k
, int l
, double tau
, double s
)
53 a
[i
][j
] = g
- s
* (h
+ g
* tau
);
54 a
[k
][l
] = h
+ s
* (g
- h
* tau
);
57 void jacobi(double** a
, int n
, double d
[], double** v
, int* nrot
)
61 double tresh
, theta
, tau
, t
, sm
, s
, h
, g
, c
, *b
, *z
;
65 for (ip
= 0; ip
< n
; ip
++)
67 for (iq
= 0; iq
< n
; iq
++)
73 for (ip
= 0; ip
< n
; ip
++)
75 b
[ip
] = d
[ip
] = a
[ip
][ip
];
79 for (i
= 1; i
<= 50; i
++)
82 for (ip
= 0; ip
< n
- 1; ip
++)
84 for (iq
= ip
+ 1; iq
< n
; iq
++)
86 sm
+= std::abs(a
[ip
][iq
]);
97 tresh
= 0.2 * sm
/ (n
* n
);
103 for (ip
= 0; ip
< n
- 1; ip
++)
105 for (iq
= ip
+ 1; iq
< n
; iq
++)
107 g
= 100.0 * std::abs(a
[ip
][iq
]);
108 if (i
> 4 && std::abs(d
[ip
]) + g
== std::abs(d
[ip
])
109 && std::abs(d
[iq
]) + g
== std::abs(d
[iq
]))
113 else if (std::abs(a
[ip
][iq
]) > tresh
)
116 if (std::abs(h
) + g
== std::abs(h
))
122 theta
= 0.5 * h
/ (a
[ip
][iq
]);
123 t
= 1.0 / (std::abs(theta
) + std::sqrt(1.0 + theta
* theta
));
129 c
= 1.0 / std::sqrt(1 + t
* t
);
138 for (j
= 0; j
< ip
; j
++)
140 do_rotate(a
, j
, ip
, j
, iq
, tau
, s
);
142 for (j
= ip
+ 1; j
< iq
; j
++)
144 do_rotate(a
, ip
, j
, j
, iq
, tau
, s
);
146 for (j
= iq
+ 1; j
< n
; j
++)
148 do_rotate(a
, ip
, j
, iq
, j
, tau
, s
);
150 for (j
= 0; j
< n
; j
++)
152 do_rotate(v
, j
, ip
, j
, iq
, tau
, s
);
158 for (ip
= 0; ip
< n
; ip
++)
165 gmx_fatal(FARGS
, "Error: Too many iterations in routine JACOBI\n");
168 int m_inv_gen(real
* m
, int n
, real
* minv
)
170 double **md
, **v
, *eig
, tol
, s
;
171 int nzero
, i
, j
, k
, nrot
;
174 for (i
= 0; i
< n
; i
++)
179 for (i
= 0; i
< n
; i
++)
184 for (i
= 0; i
< n
; i
++)
186 for (j
= 0; j
< n
; j
++)
188 md
[i
][j
] = m
[i
* n
+ j
];
193 for (i
= 0; i
< n
; i
++)
195 tol
+= std::abs(md
[i
][i
]);
197 tol
= 1e-6 * tol
/ n
;
199 jacobi(md
, n
, eig
, v
, &nrot
);
202 for (i
= 0; i
< n
; i
++)
204 if (std::abs(eig
[i
]) < tol
)
211 eig
[i
] = 1.0 / eig
[i
];
215 for (i
= 0; i
< n
; i
++)
217 for (j
= 0; j
< n
; j
++)
220 for (k
= 0; k
< n
; k
++)
222 s
+= eig
[k
] * v
[i
][k
] * v
[j
][k
];
229 for (i
= 0; i
< n
; i
++)
234 for (i
= 0; i
< n
; i
++)