Split lines with many copyright years
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / xlate.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-2004, 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 "xlate.h"
41
42 #include <cctype>
43 #include <cstring>
44
45 #include <string>
46 #include <vector>
47
48 #include "gromacs/gmxpreprocess/fflibutil.h"
49 #include "gromacs/gmxpreprocess/grompp_impl.h"
50 #include "gromacs/topology/residuetypes.h"
51 #include "gromacs/topology/symtab.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
56 #include "gromacs/utility/strdb.h"
57
58 #include "hackblock.h"
59
60 typedef struct
61 {
62     char* filebase;
63     char* res;
64     char* atom;
65     char* replace;
66 } t_xlate_atom;
67
68 static void get_xlatoms(const std::string& filename, FILE* fp, int* nptr, t_xlate_atom** xlptr)
69 {
70     char          filebase[STRLEN];
71     char          line[STRLEN];
72     char          abuf[1024], rbuf[1024], repbuf[1024], dumbuf[1024];
73     char*         _ptr;
74     int           n, na, idum;
75     t_xlate_atom* xl;
76
77     fflib_filename_base(filename.c_str(), filebase, STRLEN);
78
79     n  = *nptr;
80     xl = *xlptr;
81
82     while (get_a_line(fp, line, STRLEN))
83     {
84         na = sscanf(line, "%s%s%s%s", rbuf, abuf, repbuf, dumbuf);
85         /* Check if we are reading an old format file with the number of items
86          * on the first line.
87          */
88         if (na == 1 && n == *nptr && sscanf(rbuf, "%d", &idum) == 1)
89         {
90             continue;
91         }
92         if (na != 3)
93         {
94             gmx_fatal(FARGS, "Expected a residue name and two atom names in file '%s', not '%s'",
95                       filename.c_str(), line);
96         }
97
98         srenew(xl, n + 1);
99         xl[n].filebase = gmx_strdup(filebase);
100
101         /* Use wildcards... */
102         if (strcmp(rbuf, "*") != 0)
103         {
104             xl[n].res = gmx_strdup(rbuf);
105         }
106         else
107         {
108             xl[n].res = nullptr;
109         }
110
111         /* Replace underscores in the string by spaces */
112         while ((_ptr = strchr(abuf, '_')) != nullptr)
113         {
114             *_ptr = ' ';
115         }
116
117         xl[n].atom    = gmx_strdup(abuf);
118         xl[n].replace = gmx_strdup(repbuf);
119         n++;
120     }
121
122     *nptr  = n;
123     *xlptr = xl;
124 }
125
126 static void done_xlatom(int nxlate, t_xlate_atom* xlatom)
127 {
128     int i;
129
130     for (i = 0; (i < nxlate); i++)
131     {
132         sfree(xlatom[i].filebase);
133         if (xlatom[i].res != nullptr)
134         {
135             sfree(xlatom[i].res);
136         }
137         sfree(xlatom[i].atom);
138         sfree(xlatom[i].replace);
139     }
140     sfree(xlatom);
141 }
142
143 void rename_atoms(const char*                            xlfile,
144                   const char*                            ffdir,
145                   t_atoms*                               atoms,
146                   t_symtab*                              symtab,
147                   gmx::ArrayRef<const PreprocessResidue> localPpResidue,
148                   bool                                   bResname,
149                   ResidueType*                           rt,
150                   bool                                   bReorderNum,
151                   bool                                   bVerbose)
152 {
153     int           nxlate, a, i, resind;
154     t_xlate_atom* xlatom;
155     char          c, *rnm, atombuf[32];
156     bool          bReorderedNum, bRenamed, bMatch;
157     bool          bStartTerm, bEndTerm;
158
159     nxlate = 0;
160     xlatom = nullptr;
161     if (xlfile != nullptr)
162     {
163         gmx::FilePtr fp = gmx::openLibraryFile(xlfile);
164         get_xlatoms(xlfile, fp.get(), &nxlate, &xlatom);
165     }
166     else
167     {
168         std::vector<std::string> fns = fflib_search_file_end(ffdir, ".arn", FALSE);
169         for (const auto& filename : fns)
170         {
171             FILE* fp = fflib_open(filename);
172             get_xlatoms(filename, fp, &nxlate, &xlatom);
173             gmx_ffclose(fp);
174         }
175     }
176
177     for (a = 0; (a < atoms->nr); a++)
178     {
179         resind = atoms->atom[a].resind;
180
181         bStartTerm = (resind == 0) || atoms->resinfo[resind].chainnum != atoms->resinfo[resind - 1].chainnum;
182         bEndTerm = (resind >= atoms->nres - 1)
183                    || atoms->resinfo[resind].chainnum != atoms->resinfo[resind + 1].chainnum;
184
185         if (bResname)
186         {
187             rnm = *(atoms->resinfo[resind].name);
188         }
189         else
190         {
191             rnm = *(atoms->resinfo[resind].rtp);
192         }
193
194         strcpy(atombuf, *(atoms->atomname[a]));
195         bReorderedNum = FALSE;
196         if (bReorderNum)
197         {
198             if (isdigit(atombuf[0]))
199             {
200                 c = atombuf[0];
201                 for (i = 0; (static_cast<size_t>(i) < strlen(atombuf) - 1); i++)
202                 {
203                     atombuf[i] = atombuf[i + 1];
204                 }
205                 atombuf[i]    = c;
206                 bReorderedNum = TRUE;
207             }
208         }
209         bRenamed = FALSE;
210         for (i = 0; (i < nxlate) && !bRenamed; i++)
211         {
212             /* Check if the base file name of the rtp and arn entry match */
213             if (localPpResidue.empty()
214                 || gmx::equalCaseInsensitive(localPpResidue[resind].filebase, xlatom[i].filebase))
215             {
216                 /* Match the residue name */
217                 bMatch = (xlatom[i].res == nullptr
218                           || (gmx_strcasecmp("protein-nterm", xlatom[i].res) == 0
219                               && rt->namedResidueHasType(rnm, "Protein") && bStartTerm)
220                           || (gmx_strcasecmp("protein-cterm", xlatom[i].res) == 0
221                               && rt->namedResidueHasType(rnm, "Protein") && bEndTerm)
222                           || (gmx_strcasecmp("protein", xlatom[i].res) == 0
223                               && rt->namedResidueHasType(rnm, "Protein"))
224                           || (gmx_strcasecmp("DNA", xlatom[i].res) == 0
225                               && rt->namedResidueHasType(rnm, "DNA"))
226                           || (gmx_strcasecmp("RNA", xlatom[i].res) == 0
227                               && rt->namedResidueHasType(rnm, "RNA")));
228                 if (!bMatch)
229                 {
230                     const char* ptr0 = rnm;
231                     const char* ptr1 = xlatom[i].res;
232                     while (ptr0[0] != '\0' && ptr1[0] != '\0' && (ptr0[0] == ptr1[0] || ptr1[0] == '?'))
233                     {
234                         ptr0++;
235                         ptr1++;
236                     }
237                     bMatch = (ptr0[0] == '\0' && ptr1[0] == '\0');
238                 }
239                 if (bMatch && strcmp(atombuf, xlatom[i].atom) == 0)
240                 {
241                     /* We have a match. */
242                     /* Don't free the old atomname,
243                      * since it might be in the symtab.
244                      */
245                     const char* ptr0 = xlatom[i].replace;
246                     if (bVerbose)
247                     {
248                         printf("Renaming atom '%s' in residue %d %s to '%s'\n", *atoms->atomname[a],
249                                atoms->resinfo[resind].nr, *atoms->resinfo[resind].name, ptr0);
250                     }
251                     atoms->atomname[a] = put_symtab(symtab, ptr0);
252                     bRenamed           = TRUE;
253                 }
254             }
255         }
256         if (bReorderedNum && !bRenamed)
257         {
258             atoms->atomname[a] = put_symtab(symtab, atombuf);
259         }
260     }
261
262     done_xlatom(nxlate, xlatom);
263 }