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