Remove unused thole polarization rfac parameter
[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,2021, 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/mdtypes/atominfo.h"
45 #include "gromacs/topology/ifunc.h"
46 #include "gromacs/topology/topology.h"
47 #include "gromacs/utility/arrayref.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/smalloc.h"
50
51 static gmx_bool ip_pert(int ftype, const t_iparams* ip)
52 {
53     if (NRFPB(ftype) == 0)
54     {
55         return FALSE;
56     }
57
58     bool bPert = false;
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 (int 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 (int 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,
124                       "Function type %s does not support currentely free energy calculations",
125                       interaction_function[ftype].longname);
126         default:
127             gmx_fatal(FARGS,
128                       "Function type %s not implemented in ip_pert",
129                       interaction_function[ftype].longname);
130     }
131
132     return bPert;
133 }
134
135
136 //! Return whether the two atom indices in a 1-4 interaction have perturbed charges per \c atomInfo
137 static bool hasPerturbedChargeIn14Interaction(int atom0, int atom1, gmx::ArrayRef<const int64_t> atomInfo)
138 {
139     return (bool(atomInfo[atom0] & gmx::sc_atomInfo_HasPerturbedChargeIn14Interaction)
140             || bool(atomInfo[atom1] & gmx::sc_atomInfo_HasPerturbedChargeIn14Interaction));
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, gmx::ArrayRef<const int64_t> atomInfo)
180 {
181     bool havePerturbedInteractions = false;
182
183     int      iabuf_nalloc = 0;
184     t_iatom* iabuf        = nullptr;
185
186     for (int ftype = 0; ftype < F_NRE; ftype++)
187     {
188         if (interaction_function[ftype].flags & IF_BOND)
189         {
190             InteractionList* ilist  = &idef->il[ftype];
191             int*             iatoms = ilist->iatoms.data();
192             const int        nral   = NRAL(ftype);
193             int              ic     = 0;
194             int              ib     = 0;
195             int              i      = 0;
196             while (i < ilist->size())
197             {
198                 /* Check if this interaction is perturbed */
199                 if (ip_pert(ftype, idef->iparams.data() + iatoms[i])
200                     || (ftype == F_LJ14
201                         && hasPerturbedChargeIn14Interaction(iatoms[i + 1], iatoms[i + 2], atomInfo)))
202                 {
203                     /* Copy to the perturbed buffer */
204                     if (ib + 1 + nral > iabuf_nalloc)
205                     {
206                         iabuf_nalloc = over_alloc_large(ib + 1 + nral);
207                         srenew(iabuf, iabuf_nalloc);
208                     }
209                     for (int a = 0; a < 1 + nral; a++)
210                     {
211                         iabuf[ib++] = iatoms[i++];
212                     }
213
214                     havePerturbedInteractions = true;
215                 }
216                 else
217                 {
218                     /* Copy in place */
219                     for (int a = 0; a < 1 + nral; a++)
220                     {
221                         iatoms[ic++] = iatoms[i++];
222                     }
223                 }
224             }
225             /* Now we know the number of non-perturbed interactions */
226             idef->numNonperturbedInteractions[ftype] = ic;
227
228             /* Copy the buffer with perturbed interactions to the ilist */
229             for (int a = 0; a < ib; a++)
230             {
231                 iatoms[ic++] = iabuf[a];
232             }
233
234             if (debug)
235             {
236                 const int numNonperturbed = idef->numNonperturbedInteractions[ftype];
237                 fprintf(debug,
238                         "%s non-pert %d pert %d\n",
239                         interaction_function[ftype].longname,
240                         numNonperturbed,
241                         ilist->size() - numNonperturbed);
242             }
243         }
244     }
245
246     sfree(iabuf);
247
248     idef->ilsort = (havePerturbedInteractions ? ilsortFE_SORTED : ilsortNO_FE);
249 }