Convert gmx_ffparams_t to C++
[alexxy/gromacs.git] / src / gromacs / topology / topsort.cpp
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-2008, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2017,2018, 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 "topsort.h"
40
41 #include <cstdio>
42
43 #include "gromacs/topology/ifunc.h"
44 #include "gromacs/topology/topology.h"
45 #include "gromacs/utility/arrayref.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/smalloc.h"
48
49 static gmx_bool ip_pert(int ftype, const t_iparams *ip)
50 {
51     gmx_bool bPert;
52     int      i;
53
54     if (NRFPB(ftype) == 0)
55     {
56         return FALSE;
57     }
58
59     switch (ftype)
60     {
61         case F_BONDS:
62         case F_G96BONDS:
63         case F_HARMONIC:
64         case F_ANGLES:
65         case F_G96ANGLES:
66         case F_IDIHS:
67             bPert = (ip->harmonic.rA  != ip->harmonic.rB ||
68                      ip->harmonic.krA != ip->harmonic.krB);
69             break;
70         case F_MORSE:
71             bPert = (ip->morse.b0A  != ip->morse.b0B ||
72                      ip->morse.cbA  != ip->morse.cbB ||
73                      ip->morse.betaA  != ip->morse.betaB);
74             break;
75         case F_RESTRBONDS:
76             bPert = (ip->restraint.lowA  != ip->restraint.lowB ||
77                      ip->restraint.up1A  != ip->restraint.up1B ||
78                      ip->restraint.up2A  != ip->restraint.up2B ||
79                      ip->restraint.kA    != ip->restraint.kB);
80             break;
81         case F_UREY_BRADLEY:
82             bPert = (ip->u_b.thetaA  != ip->u_b.thetaB  ||
83                      ip->u_b.kthetaA != ip->u_b.kthetaB ||
84                      ip->u_b.r13A    != ip->u_b.r13B    ||
85                      ip->u_b.kUBA    != ip->u_b.kUBB);
86             break;
87         case F_PDIHS:
88         case F_PIDIHS:
89         case F_ANGRES:
90         case F_ANGRESZ:
91             bPert = (ip->pdihs.phiA != ip->pdihs.phiB ||
92                      ip->pdihs.cpA  != ip->pdihs.cpB);
93             break;
94         case F_RBDIHS:
95             bPert = FALSE;
96             for (i = 0; i < NR_RBDIHS; i++)
97             {
98                 if (ip->rbdihs.rbcA[i] != ip->rbdihs.rbcB[i])
99                 {
100                     bPert = TRUE;
101                 }
102             }
103             break;
104         case F_TABBONDS:
105         case F_TABBONDSNC:
106         case F_TABANGLES:
107         case F_TABDIHS:
108             bPert = (ip->tab.kA != ip->tab.kB);
109             break;
110         case F_POSRES:
111             bPert = FALSE;
112             for (i = 0; i < DIM; i++)
113             {
114                 if (ip->posres.pos0A[i] != ip->posres.pos0B[i] ||
115                     ip->posres.fcA[i]   != ip->posres.fcB[i])
116                 {
117                     bPert = TRUE;
118                 }
119             }
120             break;
121         case F_DIHRES:
122             bPert = ((ip->dihres.phiA != ip->dihres.phiB) ||
123                      (ip->dihres.dphiA != ip->dihres.dphiB) ||
124                      (ip->dihres.kfacA != ip->dihres.kfacB));
125             break;
126         case F_LJ14:
127             bPert = (ip->lj14.c6A  != ip->lj14.c6B ||
128                      ip->lj14.c12A != ip->lj14.c12B);
129             break;
130         case F_CMAP:
131             bPert = FALSE;
132             break;
133         case F_RESTRANGLES:
134         case F_RESTRDIHS:
135         case F_CBTDIHS:
136             gmx_fatal(FARGS, "Function type %s does not support currentely free energy calculations",
137                       interaction_function[ftype].longname);
138         default:
139             gmx_fatal(FARGS, "Function type %s not implemented in ip_pert",
140                       interaction_function[ftype].longname);
141     }
142
143     return bPert;
144 }
145
146 static gmx_bool ip_q_pert(int ftype, const t_iatom *ia,
147                           const t_iparams *ip, const real *qA, const real *qB)
148 {
149     /* 1-4 interactions do not have the charges stored in the iparams list,
150      * so we need a separate check for those.
151      */
152     return (ip_pert(ftype, ip+ia[0]) ||
153             (ftype == F_LJ14 && (qA[ia[1]] != qB[ia[1]] ||
154                                  qA[ia[2]] != qB[ia[2]])));
155 }
156
157 gmx_bool gmx_mtop_bondeds_free_energy(const gmx_mtop_t *mtop)
158 {
159     const gmx_ffparams_t *ffparams = &mtop->ffparams;
160
161     /* Loop over all the function types and compare the A/B parameters */
162     gmx_bool bPert = FALSE;
163     for (int i = 0; i < ffparams->numTypes(); i++)
164     {
165         int ftype = ffparams->functype[i];
166         if (interaction_function[ftype].flags & IF_BOND)
167         {
168             if (ip_pert(ftype, &ffparams->iparams[i]))
169             {
170                 bPert = TRUE;
171             }
172         }
173     }
174
175     /* Check perturbed charges for 1-4 interactions */
176     for (const gmx_molblock_t &molb : mtop->molblock)
177     {
178         const t_atom             *atom = mtop->moltype[molb.type].atoms.atom;
179         const InteractionList    &il   = mtop->moltype[molb.type].ilist[F_LJ14];
180         gmx::ArrayRef<const int>  ia   = il.iatoms;
181         for (int i = 0; i < il.size(); i += 3)
182         {
183             if (atom[ia[i+1]].q != atom[ia[i+1]].qB ||
184                 atom[ia[i+2]].q != atom[ia[i+2]].qB)
185             {
186                 bPert = TRUE;
187             }
188         }
189     }
190
191     return bPert;
192 }
193
194 void gmx_sort_ilist_fe(t_idef *idef, const real *qA, const real *qB)
195 {
196     int        ftype, nral, i, ic, ib, a;
197     t_ilist   *ilist;
198     t_iatom   *iatoms;
199     t_iatom   *iabuf;
200     int        iabuf_nalloc;
201
202     if (qB == nullptr)
203     {
204         qB = qA;
205     }
206
207     iabuf_nalloc = 0;
208     iabuf        = nullptr;
209
210     const t_iparams *iparams = idef->iparams;
211
212     for (ftype = 0; ftype < F_NRE; ftype++)
213     {
214         if (interaction_function[ftype].flags & IF_BOND)
215         {
216             ilist  = &idef->il[ftype];
217             iatoms = ilist->iatoms;
218             nral   = NRAL(ftype);
219             ic     = 0;
220             ib     = 0;
221             i      = 0;
222             while (i < ilist->nr)
223             {
224                 /* Check if this interaction is perturbed */
225                 if (ip_q_pert(ftype, iatoms+i, iparams, qA, qB))
226                 {
227                     /* Copy to the perturbed buffer */
228                     if (ib + 1 + nral > iabuf_nalloc)
229                     {
230                         iabuf_nalloc = over_alloc_large(ib+1+nral);
231                         srenew(iabuf, iabuf_nalloc);
232                     }
233                     for (a = 0; a < 1+nral; a++)
234                     {
235                         iabuf[ib++] = iatoms[i++];
236                     }
237                 }
238                 else
239                 {
240                     /* Copy in place */
241                     for (a = 0; a < 1+nral; a++)
242                     {
243                         iatoms[ic++] = iatoms[i++];
244                     }
245                 }
246             }
247             /* Now we now the number of non-perturbed interactions */
248             ilist->nr_nonperturbed = ic;
249
250             /* Copy the buffer with perturbed interactions to the ilist */
251             for (a = 0; a < ib; a++)
252             {
253                 iatoms[ic++] = iabuf[a];
254             }
255
256             if (debug)
257             {
258                 fprintf(debug, "%s non-pert %d pert %d\n",
259                         interaction_function[ftype].longname,
260                         ilist->nr_nonperturbed,
261                         ilist->nr-ilist->nr_nonperturbed);
262             }
263         }
264     }
265
266     sfree(iabuf);
267
268     idef->ilsort = ilsortFE_SORTED;
269 }