Split lines with many copyright years
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / resall.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,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,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 "resall.h"
41
42 #include <cctype>
43 #include <cstdlib>
44 #include <cstring>
45
46 #include <algorithm>
47 #include <string>
48 #include <vector>
49
50 #include "gromacs/gmxpreprocess/fflibutil.h"
51 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
52 #include "gromacs/gmxpreprocess/grompp_impl.h"
53 #include "gromacs/gmxpreprocess/notset.h"
54 #include "gromacs/gmxpreprocess/pgutil.h"
55 #include "gromacs/topology/atoms.h"
56 #include "gromacs/topology/symtab.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
61 #include "gromacs/utility/strdb.h"
62
63 #include "hackblock.h"
64
65 PreprocessingAtomTypes read_atype(const char* ffdir, t_symtab* tab)
66 {
67     FILE*   in;
68     char    buf[STRLEN], name[STRLEN];
69     double  m;
70     t_atom* a;
71
72     std::vector<std::string> files = fflib_search_file_end(ffdir, ".atp", TRUE);
73     snew(a, 1);
74     PreprocessingAtomTypes at;
75
76     for (const auto& filename : files)
77     {
78         in = fflib_open(filename);
79         while (!feof(in))
80         {
81             /* Skip blank or comment-only lines */
82             do
83             {
84                 if (fgets2(buf, STRLEN, in) != nullptr)
85                 {
86                     strip_comment(buf);
87                     trim(buf);
88                 }
89             } while ((feof(in) == 0) && strlen(buf) == 0);
90
91             if (sscanf(buf, "%s%lf", name, &m) == 2)
92             {
93                 a->m = m;
94                 at.addType(tab, *a, name, InteractionOfType({}, {}), 0, 0);
95                 fflush(stderr);
96             }
97             else
98             {
99                 gmx_fatal(FARGS, "Invalid atomtype format: '%s'", buf);
100             }
101         }
102         gmx_ffclose(in);
103     }
104     sfree(a);
105     return at;
106 }
107
108 static void print_resatoms(FILE* out, const PreprocessingAtomTypes& atype, const PreprocessResidue& rtpDBEntry)
109 {
110     fprintf(out, "[ %s ]\n", rtpDBEntry.resname.c_str());
111     fprintf(out, " [ atoms ]\n");
112
113     for (int j = 0; (j < rtpDBEntry.natom()); j++)
114     {
115         int         tp   = rtpDBEntry.atom[j].type;
116         const char* tpnm = atype.atomNameFromAtomType(tp);
117         if (tpnm == nullptr)
118         {
119             gmx_fatal(FARGS, "Incorrect atomtype (%d)", tp);
120         }
121         fprintf(out, "%6s  %6s  %8.3f  %6d\n", *(rtpDBEntry.atomname[j]), tpnm,
122                 rtpDBEntry.atom[j].q, rtpDBEntry.cgnr[j]);
123     }
124 }
125
126 static bool read_atoms(FILE* in, char* line, PreprocessResidue* r0, t_symtab* tab, PreprocessingAtomTypes* atype)
127 {
128     int    cg;
129     char   buf[256], buf1[256];
130     double q;
131
132     /* Read Atoms */
133     r0->atom.clear();
134     r0->atomname.clear();
135     r0->cgnr.clear();
136     while (get_a_line(in, line, STRLEN) && (strchr(line, '[') == nullptr))
137     {
138         if (sscanf(line, "%s%s%lf%d", buf, buf1, &q, &cg) != 4)
139         {
140             return FALSE;
141         }
142         r0->atomname.push_back(put_symtab(tab, buf));
143         r0->atom.emplace_back();
144         r0->atom.back().q = q;
145         r0->cgnr.push_back(cg);
146         int j = atype->atomTypeFromName(buf1);
147         if (j == NOTSET)
148         {
149             gmx_fatal(FARGS,
150                       "Atom type %s (residue %s) not found in atomtype "
151                       "database",
152                       buf1, r0->resname.c_str());
153         }
154         r0->atom.back().type = j;
155         r0->atom.back().m    = atype->atomMassFromAtomType(j);
156     }
157
158     return TRUE;
159 }
160
161 static bool read_bondeds(int bt, FILE* in, char* line, PreprocessResidue* rtpDBEntry)
162 {
163     char str[STRLEN];
164
165     while (get_a_line(in, line, STRLEN) && (strchr(line, '[') == nullptr))
166     {
167         int n = 0;
168         int ni;
169         rtpDBEntry->rb[bt].b.emplace_back();
170         BondedInteraction* newBond = &rtpDBEntry->rb[bt].b.back();
171         for (int j = 0; j < btsNiatoms[bt]; j++)
172         {
173             if (sscanf(line + n, "%s%n", str, &ni) == 1)
174             {
175                 newBond->a[j] = str;
176             }
177             else
178             {
179                 return FALSE;
180             }
181             n += ni;
182         }
183         while (isspace(line[n]))
184         {
185             n++;
186         }
187         rtrim(line + n);
188         newBond->s = line + n;
189     }
190
191     return TRUE;
192 }
193
194 static void print_resbondeds(FILE* out, int bt, const PreprocessResidue& rtpDBEntry)
195 {
196     if (!rtpDBEntry.rb[bt].b.empty())
197     {
198         fprintf(out, " [ %s ]\n", btsNames[bt]);
199
200         for (const auto& b : rtpDBEntry.rb[bt].b)
201         {
202             for (int j = 0; j < btsNiatoms[bt]; j++)
203             {
204                 fprintf(out, "%6s ", b.a[j].c_str());
205             }
206             if (!b.s.empty())
207             {
208                 fprintf(out, "    %s", b.s.c_str());
209             }
210             fprintf(out, "\n");
211         }
212     }
213 }
214
215 static void check_rtp(gmx::ArrayRef<const PreprocessResidue> rtpDBEntry, const std::string& libfn)
216 {
217     /* check for double entries, assuming list is already sorted */
218     for (auto it = rtpDBEntry.begin() + 1; it != rtpDBEntry.end(); it++)
219     {
220         auto prev = it - 1;
221         if (gmx::equalCaseInsensitive(prev->resname, it->resname))
222         {
223             fprintf(stderr, "WARNING double entry %s in file %s\n", it->resname.c_str(), libfn.c_str());
224         }
225     }
226 }
227
228 static int get_bt(char* header)
229 {
230     int i;
231
232     for (i = 0; i < ebtsNR; i++)
233     {
234         if (gmx_strcasecmp(btsNames[i], header) == 0)
235         {
236             return i;
237         }
238     }
239     return NOTSET;
240 }
241
242 /* print all the ebtsNR type numbers */
243 static void print_resall_header(FILE* out, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
244 {
245     fprintf(out, "[ bondedtypes ]\n");
246     fprintf(out,
247             "; bonds  angles  dihedrals  impropers all_dihedrals nr_exclusions  HH14  "
248             "remove_dih\n");
249     fprintf(out, " %5d  %6d  %9d  %9d  %14d  %14d %14d %14d\n\n", rtpDBEntry[0].rb[0].type,
250             rtpDBEntry[0].rb[1].type, rtpDBEntry[0].rb[2].type, rtpDBEntry[0].rb[3].type,
251             static_cast<int>(rtpDBEntry[0].bKeepAllGeneratedDihedrals), rtpDBEntry[0].nrexcl,
252             static_cast<int>(rtpDBEntry[0].bGenerateHH14Interactions),
253             static_cast<int>(rtpDBEntry[0].bRemoveDihedralIfWithImproper));
254 }
255
256 void print_resall(FILE* out, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry, const PreprocessingAtomTypes& atype)
257 {
258     if (rtpDBEntry.empty())
259     {
260         return;
261     }
262
263     print_resall_header(out, rtpDBEntry);
264
265     for (const auto& r : rtpDBEntry)
266     {
267         if (r.natom() > 0)
268         {
269             print_resatoms(out, atype, r);
270             for (int bt = 0; bt < ebtsNR; bt++)
271             {
272                 print_resbondeds(out, bt, r);
273             }
274         }
275     }
276 }
277
278 void readResidueDatabase(const std::string&              rrdb,
279                          std::vector<PreprocessResidue>* rtpDBEntry,
280                          PreprocessingAtomTypes*         atype,
281                          t_symtab*                       tab,
282                          bool                            bAllowOverrideRTP)
283 {
284     FILE* in;
285     char  filebase[STRLEN], line[STRLEN], header[STRLEN];
286     int   bt, nparam;
287     int   dum1, dum2, dum3;
288     bool  bNextResidue, bError;
289
290     fflib_filename_base(rrdb.c_str(), filebase, STRLEN);
291
292     in = fflib_open(rrdb);
293
294     PreprocessResidue header_settings;
295
296     /* these bonded parameters will overwritten be when  *
297      * there is a [ bondedtypes ] entry in the .rtp file */
298     header_settings.rb[ebtsBONDS].type  = 1; /* normal bonds     */
299     header_settings.rb[ebtsANGLES].type = 1; /* normal angles    */
300     header_settings.rb[ebtsPDIHS].type  = 1; /* normal dihedrals */
301     header_settings.rb[ebtsIDIHS].type  = 2; /* normal impropers */
302     header_settings.rb[ebtsEXCLS].type  = 1; /* normal exclusions */
303     header_settings.rb[ebtsCMAP].type   = 1; /* normal cmap torsions */
304
305     header_settings.bKeepAllGeneratedDihedrals    = FALSE;
306     header_settings.nrexcl                        = 3;
307     header_settings.bGenerateHH14Interactions     = TRUE;
308     header_settings.bRemoveDihedralIfWithImproper = TRUE;
309
310     /* Column 5 & 6 aren't really bonded types, but we include
311      * them here to avoid introducing a new section:
312      * Column 5 : This controls the generation of dihedrals from the bonding.
313      *            All possible dihedrals are generated automatically. A value of
314      *            1 here means that all these are retained. A value of
315      *            0 here requires generated dihedrals be removed if
316      *              * there are any dihedrals on the same central atoms
317      *                specified in the residue topology, or
318      *              * there are other identical generated dihedrals
319      *                sharing the same central atoms, or
320      *              * there are other generated dihedrals sharing the
321      *                same central bond that have fewer hydrogen atoms
322      * Column 6: Number of bonded neighbors to exclude.
323      * Column 7: Generate 1,4 interactions between two hydrogen atoms
324      * Column 8: Remove proper dihedrals if centered on the same bond
325      *           as an improper dihedral
326      */
327     get_a_line(in, line, STRLEN);
328     if (!get_header(line, header))
329     {
330         gmx_fatal(FARGS, "in .rtp file at line:\n%s\n", line);
331     }
332     if (gmx::equalCaseInsensitive("bondedtypes", header, 5))
333     {
334         get_a_line(in, line, STRLEN);
335         if ((nparam = sscanf(line, "%d %d %d %d %d %d %d %d", &header_settings.rb[ebtsBONDS].type,
336                              &header_settings.rb[ebtsANGLES].type,
337                              &header_settings.rb[ebtsPDIHS].type, &header_settings.rb[ebtsIDIHS].type,
338                              &dum1, &header_settings.nrexcl, &dum2, &dum3))
339             < 4)
340         {
341             gmx_fatal(FARGS, "need 4 to 8 parameters in the header of .rtp file %s at line:\n%s\n",
342                       rrdb.c_str(), line);
343         }
344         header_settings.bKeepAllGeneratedDihedrals    = (dum1 != 0);
345         header_settings.bGenerateHH14Interactions     = (dum2 != 0);
346         header_settings.bRemoveDihedralIfWithImproper = (dum3 != 0);
347         get_a_line(in, line, STRLEN);
348         if (nparam < 5)
349         {
350             fprintf(stderr, "Using default: not generating all possible dihedrals\n");
351             header_settings.bKeepAllGeneratedDihedrals = FALSE;
352         }
353         if (nparam < 6)
354         {
355             fprintf(stderr, "Using default: excluding 3 bonded neighbors\n");
356             header_settings.nrexcl = 3;
357         }
358         if (nparam < 7)
359         {
360             fprintf(stderr, "Using default: generating 1,4 H--H interactions\n");
361             header_settings.bGenerateHH14Interactions = TRUE;
362         }
363         if (nparam < 8)
364         {
365             fprintf(stderr,
366                     "Using default: removing proper dihedrals found on the same bond as a proper "
367                     "dihedral\n");
368             header_settings.bRemoveDihedralIfWithImproper = TRUE;
369         }
370     }
371     else
372     {
373         fprintf(stderr,
374                 "Reading .rtp file without '[ bondedtypes ]' directive,\n"
375                 "Will proceed as if the entry was:\n");
376         print_resall_header(stderr, gmx::arrayRefFromArray(&header_settings, 1));
377     }
378     /* We don't know the current size of rrtp, but simply realloc immediately */
379     auto oldArrayEnd = rtpDBEntry->end();
380     while (!feof(in))
381     {
382         /* Initialise rtp entry structure */
383         rtpDBEntry->push_back(header_settings);
384         PreprocessResidue* res = &rtpDBEntry->back();
385         if (!get_header(line, header))
386         {
387             gmx_fatal(FARGS, "in .rtp file at line:\n%s\n", line);
388         }
389         res->resname  = header;
390         res->filebase = filebase;
391
392         get_a_line(in, line, STRLEN);
393         bError       = FALSE;
394         bNextResidue = FALSE;
395         do
396         {
397             if (!get_header(line, header))
398             {
399                 bError = TRUE;
400             }
401             else
402             {
403                 bt = get_bt(header);
404                 if (bt != NOTSET)
405                 {
406                     /* header is an bonded directive */
407                     bError = !read_bondeds(bt, in, line, res);
408                 }
409                 else if (gmx::equalCaseInsensitive("atoms", header, 5))
410                 {
411                     /* header is the atoms directive */
412                     bError = !read_atoms(in, line, res, tab, atype);
413                 }
414                 else
415                 {
416                     /* else header must be a residue name */
417                     bNextResidue = TRUE;
418                 }
419             }
420             if (bError)
421             {
422                 gmx_fatal(FARGS, "in .rtp file in residue %s at line:\n%s\n", res->resname.c_str(), line);
423             }
424         } while ((feof(in) == 0) && !bNextResidue);
425
426         if (res->natom() == 0)
427         {
428             gmx_fatal(FARGS, "No atoms found in .rtp file in residue %s\n", res->resname.c_str());
429         }
430
431         auto found = std::find_if(rtpDBEntry->begin(), rtpDBEntry->end() - 1,
432                                   [&res](const PreprocessResidue& entry) {
433                                       return gmx::equalCaseInsensitive(entry.resname, res->resname);
434                                   });
435
436         if (found != rtpDBEntry->end() - 1)
437         {
438             if (found >= oldArrayEnd)
439             {
440                 gmx_fatal(FARGS, "Found a second entry for '%s' in '%s'", res->resname.c_str(),
441                           rrdb.c_str());
442             }
443             if (bAllowOverrideRTP)
444             {
445                 fprintf(stderr,
446                         "Found another rtp entry for '%s' in '%s', ignoring this entry and keeping "
447                         "the one from '%s.rtp'\n",
448                         res->resname.c_str(), rrdb.c_str(), found->filebase.c_str());
449                 /* We should free all the data for this entry.
450                  * The current code gives a lot of dangling pointers.
451                  */
452                 rtpDBEntry->erase(rtpDBEntry->end() - 1);
453             }
454             else
455             {
456                 gmx_fatal(FARGS,
457                           "Found rtp entries for '%s' in both '%s' and '%s'. If you want the first "
458                           "definition to override the second one, set the -rtpo option of pdb2gmx.",
459                           res->resname.c_str(), found->filebase.c_str(), rrdb.c_str());
460             }
461         }
462     }
463     gmx_ffclose(in);
464
465     std::sort(rtpDBEntry->begin(), rtpDBEntry->end(), [](const PreprocessResidue& a, const PreprocessResidue& b) {
466         return std::lexicographical_compare(
467                 a.resname.begin(), a.resname.end(), b.resname.begin(), b.resname.end(),
468                 [](const char& c1, const char& c2) { return std::toupper(c1) < std::toupper(c2); });
469     });
470
471     check_rtp(*rtpDBEntry, rrdb);
472 }
473
474 /************************************************************
475  *
476  *                  SEARCH   ROUTINES
477  *
478  ***********************************************************/
479 static bool is_sign(char c)
480 {
481     return (c == '+' || c == '-');
482 }
483
484 /* Compares if the strins match except for a sign at the end
485  * of (only) one of the two.
486  */
487 static int neq_str_sign(const char* a1, const char* a2)
488 {
489     int l1, l2, lm;
490
491     l1 = static_cast<int>(strlen(a1));
492     l2 = static_cast<int>(strlen(a2));
493     lm = std::min(l1, l2);
494
495     if (lm >= 1 && ((l1 == l2 + 1 && is_sign(a1[l1 - 1])) || (l2 == l1 + 1 && is_sign(a2[l2 - 1])))
496         && gmx::equalCaseInsensitive(a1, a2, lm))
497     {
498         return lm;
499     }
500     else
501     {
502         return 0;
503     }
504 }
505
506 std::string searchResidueDatabase(const std::string& key, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
507 {
508     int         nbest, best, besti;
509     std::string bestbuf;
510
511     nbest = 0;
512     besti = -1;
513     /* We want to match at least one character */
514     best = 1;
515     for (auto it = rtpDBEntry.begin(); it != rtpDBEntry.end(); it++)
516     {
517         if (gmx::equalCaseInsensitive(key, it->resname))
518         {
519             besti = std::distance(rtpDBEntry.begin(), it);
520             nbest = 1;
521             break;
522         }
523         else
524         {
525             /* Allow a mismatch of at most a sign character (with warning) */
526             int n = neq_str_sign(key.c_str(), it->resname.c_str());
527             if (n >= best && n + 1 >= gmx::index(key.length()) && n + 1 >= gmx::index(it->resname.length()))
528             {
529                 if (n == best)
530                 {
531                     if (nbest == 1)
532                     {
533                         bestbuf = rtpDBEntry[besti].resname;
534                     }
535                     if (nbest >= 1)
536                     {
537                         bestbuf.append(" or ");
538                         bestbuf.append(it->resname);
539                     }
540                 }
541                 else
542                 {
543                     nbest = 0;
544                 }
545                 besti = std::distance(rtpDBEntry.begin(), it);
546                 best  = n;
547                 nbest++;
548             }
549         }
550     }
551     if (nbest > 1)
552     {
553         gmx_fatal(FARGS, "Residue '%s' not found in residue topology database, looks a bit like %s",
554                   key.c_str(), bestbuf.c_str());
555     }
556     else if (besti == -1)
557     {
558         gmx_fatal(FARGS, "Residue '%s' not found in residue topology database", key.c_str());
559     }
560     if (!gmx::equalCaseInsensitive(rtpDBEntry[besti].resname, key))
561     {
562         fprintf(stderr,
563                 "\nWARNING: '%s' not found in residue topology database, "
564                 "trying to use '%s'\n\n",
565                 key.c_str(), rtpDBEntry[besti].resname.c_str());
566     }
567
568     return rtpDBEntry[besti].resname;
569 }
570
571 gmx::ArrayRef<const PreprocessResidue>::const_iterator
572 getDatabaseEntry(const std::string& rtpname, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
573 {
574     auto found = std::find_if(rtpDBEntry.begin(), rtpDBEntry.end(),
575                               [&rtpname](const PreprocessResidue& entry) {
576                                   return gmx::equalCaseInsensitive(rtpname, entry.resname);
577                               });
578     if (found == rtpDBEntry.end())
579     {
580         /* This should never happen, since searchResidueDatabase should have been called
581          * before calling getDatabaseEntry.
582          */
583         gmx_fatal(FARGS, "Residue type '%s' not found in residue topology database", rtpname.c_str());
584     }
585
586     return found;
587 }