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