2 * This file is part of the GROMACS molecular simulation package.
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,2016,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.
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.
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.
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.
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.
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.
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/gmxfio.h"
47 #include "gromacs/gmxpreprocess/gen_ad.h"
48 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
49 #include "gromacs/gmxpreprocess/gpp_nextnb.h"
50 #include "gromacs/gmxpreprocess/grompp_impl.h"
51 #include "gromacs/gmxpreprocess/nm2type.h"
52 #include "gromacs/gmxpreprocess/notset.h"
53 #include "gromacs/gmxpreprocess/pdb2top.h"
54 #include "gromacs/gmxpreprocess/toppush.h"
55 #include "gromacs/gmxpreprocess/toputil.h"
56 #include "gromacs/listed_forces/bonded.h"
57 #include "gromacs/math/units.h"
58 #include "gromacs/math/utilities.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/math/vecdump.h"
61 #include "gromacs/mdtypes/md_enums.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/topology/symtab.h"
64 #include "gromacs/topology/topology.h"
65 #include "gromacs/utility/arraysize.h"
66 #include "gromacs/utility/cstringutil.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/gmxassert.h"
69 #include "gromacs/utility/smalloc.h"
71 #include "hackblock.h"
73 #define MARGIN_FAC 1.1
75 static bool is_bond(int nnm
, t_nm2type nmt
[], char *ai
, char *aj
, real blen
)
79 for (i
= 0; (i
< nnm
); i
++)
81 for (j
= 0; (j
< nmt
[i
].nbonds
); j
++)
83 if ((((gmx::equalCaseInsensitive(ai
, nmt
[i
].elem
, 1)) &&
84 (gmx::equalCaseInsensitive(aj
, nmt
[i
].bond
[j
], 1))) ||
85 ((gmx::equalCaseInsensitive(ai
, nmt
[i
].bond
[j
], 1)) &&
86 (gmx::equalCaseInsensitive(aj
, nmt
[i
].elem
, 1)))) &&
87 (fabs(blen
-nmt
[i
].blen
[j
]) <= 0.1*nmt
[i
].blen
[j
]))
96 static void mk_bonds(int nnm
, t_nm2type nmt
[],
97 t_atoms
*atoms
, const rvec x
[], InteractionsOfType
*bond
, int nbond
[],
98 bool bPBC
, matrix box
)
105 std::array
<real
, MAXFORCEPARAM
> forceParam
= {0.0};
108 set_pbc(&pbc
, -1, box
);
110 for (i
= 0; (i
< atoms
->nr
); i
++)
114 fprintf(stderr
, "\ratom %d", i
);
117 for (j
= i
+1; (j
< atoms
->nr
); j
++)
121 pbc_dx(&pbc
, x
[i
], x
[j
], dx
);
125 rvec_sub(x
[i
], x
[j
], dx
);
129 if (is_bond(nnm
, nmt
, *atoms
->atomname
[i
], *atoms
->atomname
[j
],
132 forceParam
[0] = std::sqrt(dx2
);
133 std::vector
<int> atoms
= {i
, j
};
134 add_param_to_list (bond
, InteractionOfType(atoms
, forceParam
));
140 fprintf(stderr
, "\ratom %d\n", i
);
144 static int *set_cgnr(t_atoms
*atoms
, bool bUsePDBcharge
, real
*qtot
, real
*mtot
)
151 snew(cgnr
, atoms
->nr
);
152 for (i
= 0; (i
< atoms
->nr
); i
++)
154 if (atoms
->pdbinfo
&& bUsePDBcharge
)
156 atoms
->atom
[i
].q
= atoms
->pdbinfo
[i
].bfac
;
158 qt
+= atoms
->atom
[i
].q
;
159 *qtot
+= atoms
->atom
[i
].q
;
160 *mtot
+= atoms
->atom
[i
].m
;
171 static void set_atom_type(PreprocessingAtomTypes
*atypes
,
174 InteractionsOfType
*bonds
,
181 snew(atoms
->atomtype
, atoms
->nr
);
182 nresolved
= nm2type(nnm
, nm2t
, tab
, atoms
, atypes
, nbonds
, bonds
);
183 if (nresolved
!= atoms
->nr
)
185 gmx_fatal(FARGS
, "Could only find a forcefield type for %d out of %d atoms",
186 nresolved
, atoms
->nr
);
189 fprintf(stderr
, "There are %zu different atom types in your sample\n",
193 static void lo_set_force_const(InteractionsOfType
*plist
, real c
[], int nrfp
, bool bRound
,
194 bool bDih
, bool bParam
)
199 for (auto ¶m
: plist
->interactionTypes
)
203 for (int j
= 0; j
< nrfp
; j
++)
212 sprintf(buf
, "%.2e", param
.c0());
213 sscanf(buf
, "%lf", &cc
);
223 c
[0] = (static_cast<int>(c
[0] + 3600)) % 360;
228 /* To put the minimum at the current angle rather than the maximum */
232 GMX_ASSERT(nrfp
<= MAXFORCEPARAM
/2, "Only 6 parameters may be used for an interaction");
233 std::array
<real
, MAXFORCEPARAM
> forceParam
;
234 for (int j
= 0; (j
< nrfp
); j
++)
236 forceParam
[j
] = c
[j
];
237 forceParam
[nrfp
+j
] = c
[j
];
239 param
= InteractionOfType(param
.atoms(), forceParam
);
243 static void set_force_const(gmx::ArrayRef
<InteractionsOfType
> plist
, real kb
, real kt
, real kp
, bool bRound
,
246 real c
[MAXFORCEPARAM
];
250 lo_set_force_const(&plist
[F_BONDS
], c
, 2, bRound
, FALSE
, bParam
);
252 lo_set_force_const(&plist
[F_ANGLES
], c
, 2, bRound
, FALSE
, bParam
);
255 lo_set_force_const(&plist
[F_PDIHS
], c
, 3, bRound
, TRUE
, bParam
);
258 static void calc_angles_dihs(InteractionsOfType
*ang
, InteractionsOfType
*dih
, const rvec x
[], bool bPBC
,
262 rvec r_ij
, r_kj
, r_kl
, m
, n
;
268 set_pbc(&pbc
, epbcXYZ
, box
);
270 for (auto &angle
: ang
->interactionTypes
)
275 real th
= RAD2DEG
*bond_angle(x
[ai
], x
[aj
], x
[ak
], bPBC
? &pbc
: nullptr,
276 r_ij
, r_kj
, &costh
, &t1
, &t2
);
277 angle
.setForceParameter(0, th
);
279 for (auto dihedral
: dih
->interactionTypes
)
281 int ai
= dihedral
.ai();
282 int aj
= dihedral
.aj();
283 int ak
= dihedral
.ak();
284 int al
= dihedral
.al();
285 real ph
= RAD2DEG
*dih_angle(x
[ai
], x
[aj
], x
[ak
], x
[al
], bPBC
? &pbc
: nullptr,
286 r_ij
, r_kj
, r_kl
, m
, n
, &t1
, &t2
, &t3
);
287 dihedral
.setForceParameter(0, ph
);
291 static void dump_hybridization(FILE *fp
, t_atoms
*atoms
, int nbonds
[])
293 for (int i
= 0; (i
< atoms
->nr
); i
++)
295 fprintf(fp
, "Atom %5s has %1d bonds\n", *atoms
->atomname
[i
], nbonds
[i
]);
299 static void print_pl(FILE *fp
, gmx::ArrayRef
<const InteractionsOfType
> plist
, int ftp
, const char *name
,
302 if (!plist
[ftp
].interactionTypes
.empty())
305 fprintf(fp
, "[ %s ]\n", name
);
306 int nrfp
= interaction_function
[ftp
].nrfpA
;
307 for (const auto ¶m
: plist
[ftp
].interactionTypes
)
309 gmx::ArrayRef
<const int> atoms
= param
.atoms();
310 gmx::ArrayRef
<const real
> forceParam
= param
.forceParam();
311 for (const auto &atom
: atoms
)
313 fprintf(fp
, " %5s", *atomname
[atom
]);
315 for (int j
= 0; (j
< nrfp
); j
++)
317 if (forceParam
[j
] != NOTSET
)
319 fprintf(fp
, " %10.3e", forceParam
[j
]);
327 static void print_rtp(const char *filenm
,
330 gmx::ArrayRef
<const InteractionsOfType
> plist
,
331 PreprocessingAtomTypes
*atypes
,
338 fp
= gmx_fio_fopen(filenm
, "w");
339 fprintf(fp
, "; %s\n", title
);
341 fprintf(fp
, "[ %s ]\n", *atoms
->resinfo
[0].name
);
343 fprintf(fp
, "[ atoms ]\n");
344 for (i
= 0; (i
< atoms
->nr
); i
++)
346 tp
= atoms
->atom
[i
].type
;
347 if ((tpnm
= atypes
->atomNameFromAtomType(tp
)) == nullptr)
349 gmx_fatal(FARGS
, "tp = %d, i = %d in print_rtp", tp
, i
);
351 fprintf(fp
, "%-8s %12s %8.4f %5d\n",
352 *atoms
->atomname
[i
], tpnm
,
353 atoms
->atom
[i
].q
, cgnr
[i
]);
355 print_pl(fp
, plist
, F_BONDS
, "bonds", atoms
->atomname
);
356 print_pl(fp
, plist
, F_ANGLES
, "angles", atoms
->atomname
);
357 print_pl(fp
, plist
, F_PDIHS
, "dihedrals", atoms
->atomname
);
358 print_pl(fp
, plist
, F_IDIHS
, "impropers", atoms
->atomname
);
363 int gmx_x2top(int argc
, char *argv
[])
365 const char *desc
[] = {
366 "[THISMODULE] generates a primitive topology from a coordinate file.",
367 "The program assumes all hydrogens are present when defining",
368 "the hybridization from the atom name and the number of bonds.",
369 "The program can also make an [REF].rtp[ref] entry, which you can then add",
370 "to the [REF].rtp[ref] database.[PAR]",
371 "When [TT]-param[tt] is set, equilibrium distances and angles",
372 "and force constants will be printed in the topology for all",
373 "interactions. The equilibrium distances and angles are taken",
374 "from the input coordinates, the force constant are set with",
375 "command line options.",
376 "The force fields somewhat supported currently are:[PAR]",
377 "G53a5 GROMOS96 53a5 Forcefield (official distribution)[PAR]",
378 "oplsaa OPLS-AA/L all-atom force field (2001 aminoacid dihedrals)[PAR]",
379 "The corresponding data files can be found in the library directory",
380 "with name [TT]atomname2type.n2t[tt]. Check Chapter 5 of the manual for more",
381 "information about file formats. By default, the force field selection",
382 "is interactive, but you can use the [TT]-ff[tt] option to specify",
383 "one of the short names above on the command line instead. In that",
384 "case [THISMODULE] just looks for the corresponding file.[PAR]",
386 const char *bugs
[] = {
387 "The atom type selection is primitive. Virtually no chemical knowledge is used",
388 "Periodic boundary conditions screw up the bonding",
389 "No improper dihedrals are generated",
390 "The atoms to atomtype translation table is incomplete ([TT]atomname2type.n2t[tt] file in the data directory). Please extend it and send the results back to the GROMACS crew."
393 std::array
<InteractionsOfType
, F_NRE
> plist
;
399 char forcefield
[32], ffdir
[STRLEN
];
400 rvec
*x
; /* coordinates? */
402 int bts
[] = { 1, 1, 1, 2 };
403 matrix box
; /* box length matrix */
404 int natoms
; /* number of atoms in one molecule */
406 bool bRTP
, bTOP
, bOPLS
;
410 gmx_output_env_t
*oenv
;
413 { efSTX
, "-f", "conf", ffREAD
},
414 { efTOP
, "-o", "out", ffOPTWR
},
415 { efRTP
, "-r", "out", ffOPTWR
}
417 #define NFILE asize(fnm)
418 real kb
= 4e5
, kt
= 400, kp
= 5;
419 PreprocessResidue rtp_header_settings
;
420 bool bRemoveDihedralIfWithImproper
= FALSE
;
421 bool bGenerateHH14Interactions
= TRUE
;
422 bool bKeepAllGeneratedDihedrals
= FALSE
;
424 bool bParam
= TRUE
, bRound
= TRUE
;
425 bool bPairs
= TRUE
, bPBC
= TRUE
;
426 bool bUsePDBcharge
= FALSE
, bVerbose
= FALSE
;
427 const char *molnm
= "ICE";
428 const char *ff
= "oplsaa";
430 { "-ff", FALSE
, etSTR
, {&ff
},
431 "Force field for your simulation. Type \"select\" for interactive selection." },
432 { "-v", FALSE
, etBOOL
, {&bVerbose
},
433 "Generate verbose output in the top file." },
434 { "-nexcl", FALSE
, etINT
, {&nrexcl
},
435 "Number of exclusions" },
436 { "-H14", FALSE
, etBOOL
, {&bGenerateHH14Interactions
},
437 "Use 3rd neighbour interactions for hydrogen atoms" },
438 { "-alldih", FALSE
, etBOOL
, {&bKeepAllGeneratedDihedrals
},
439 "Generate all proper dihedrals" },
440 { "-remdih", FALSE
, etBOOL
, {&bRemoveDihedralIfWithImproper
},
441 "Remove dihedrals on the same bond as an improper" },
442 { "-pairs", FALSE
, etBOOL
, {&bPairs
},
443 "Output 1-4 interactions (pairs) in topology file" },
444 { "-name", FALSE
, etSTR
, {&molnm
},
445 "Name of your molecule" },
446 { "-pbc", FALSE
, etBOOL
, {&bPBC
},
447 "Use periodic boundary conditions." },
448 { "-pdbq", FALSE
, etBOOL
, {&bUsePDBcharge
},
449 "Use the B-factor supplied in a [REF].pdb[ref] file for the atomic charges" },
450 { "-param", FALSE
, etBOOL
, {&bParam
},
451 "Print parameters in the output" },
452 { "-round", FALSE
, etBOOL
, {&bRound
},
453 "Round off measured values" },
454 { "-kb", FALSE
, etREAL
, {&kb
},
455 "Bonded force constant (kJ/mol/nm^2)" },
456 { "-kt", FALSE
, etREAL
, {&kt
},
457 "Angle force constant (kJ/mol/rad^2)" },
458 { "-kp", FALSE
, etREAL
, {&kp
},
459 "Dihedral angle force constant (kJ/mol/rad^2)" }
462 if (!parse_common_args(&argc
, argv
, 0, NFILE
, fnm
, asize(pa
), pa
,
463 asize(desc
), desc
, asize(bugs
), bugs
, &oenv
))
467 bRTP
= opt2bSet("-r", NFILE
, fnm
);
468 bTOP
= opt2bSet("-o", NFILE
, fnm
);
469 /* C89 requirements mean that these struct members cannot be used in
470 * the declaration of pa. So some temporary variables are needed. */
471 rtp_header_settings
.bRemoveDihedralIfWithImproper
= bRemoveDihedralIfWithImproper
;
472 rtp_header_settings
.bGenerateHH14Interactions
= bGenerateHH14Interactions
;
473 rtp_header_settings
.bKeepAllGeneratedDihedrals
= bKeepAllGeneratedDihedrals
;
474 rtp_header_settings
.nrexcl
= nrexcl
;
478 gmx_fatal(FARGS
, "Specify at least one output file");
481 /* Force field selection, interactive or direct */
482 choose_ff(strcmp(ff
, "select") == 0 ? nullptr : ff
,
483 forcefield
, sizeof(forcefield
),
484 ffdir
, sizeof(ffdir
));
486 bOPLS
= (strcmp(forcefield
, "oplsaa") == 0);
489 mymol
.name
= gmx_strdup(molnm
);
492 /* Read coordinates */
495 read_tps_conf(opt2fn("-f", NFILE
, fnm
), top
, &epbc
, &x
, nullptr, box
, FALSE
);
496 t_atoms
*atoms
= &top
->atoms
;
498 if (atoms
->pdbinfo
== nullptr)
500 snew(atoms
->pdbinfo
, natoms
);
503 sprintf(n2t
, "%s", ffdir
);
504 nm2t
= rd_nm2type(n2t
, &nnm
);
507 gmx_fatal(FARGS
, "No or incorrect atomname2type.n2t file found (looking for %s)",
512 printf("There are %d name to type translations in file %s\n", nnm
, n2t
);
516 dump_nm2type(debug
, nnm
, nm2t
);
518 printf("Generating bonds from distances...\n");
519 snew(nbonds
, atoms
->nr
);
520 mk_bonds(nnm
, nm2t
, atoms
, x
, &(plist
[F_BONDS
]), nbonds
, bPBC
, box
);
522 open_symtab(&symtab
);
523 PreprocessingAtomTypes atypes
;
524 set_atom_type(&atypes
, &symtab
, atoms
, &(plist
[F_BONDS
]), nbonds
, nnm
, nm2t
);
526 /* Make Angles and Dihedrals */
527 snew(excls
, atoms
->nr
);
528 printf("Generating angles and dihedrals from bonds...\n");
529 init_nnb(&nnb
, atoms
->nr
, 4);
530 gen_nnb(&nnb
, plist
);
531 print_nnb(&nnb
, "NNB");
532 gen_pad(&nnb
, atoms
, gmx::arrayRefFromArray(&rtp_header_settings
, 1), plist
, excls
, {}, TRUE
);
537 plist
[F_LJ14
].interactionTypes
.clear();
540 "There are %4zu %s dihedrals, %4zu impropers, %4zu angles\n"
541 " %4zu pairs, %4zu bonds and %4d atoms\n",
542 plist
[F_PDIHS
].size(),
543 bOPLS
? "Ryckaert-Bellemans" : "proper",
544 plist
[F_IDIHS
].size(), plist
[F_ANGLES
].size(),
545 plist
[F_LJ14
].size(), plist
[F_BONDS
].size(), atoms
->nr
);
547 calc_angles_dihs(&plist
[F_ANGLES
], &plist
[F_PDIHS
], x
, bPBC
, box
);
549 set_force_const(plist
, kb
, kt
, kp
, bRound
, bParam
);
551 cgnr
= set_cgnr(atoms
, bUsePDBcharge
, &qtot
, &mtot
);
552 printf("Total charge is %g, total mass is %g\n", qtot
, mtot
);
561 fp
= ftp2FILE(efTOP
, NFILE
, fnm
, "w");
562 print_top_header(fp
, ftp2fn(efTOP
, NFILE
, fnm
), TRUE
, ffdir
, 1.0);
564 write_top(fp
, nullptr, mymol
.name
.c_str(), atoms
, FALSE
, bts
, plist
, excls
, &atypes
,
565 cgnr
, rtp_header_settings
.nrexcl
);
566 print_top_mols(fp
, mymol
.name
.c_str(), ffdir
, nullptr, {}, gmx::arrayRefFromArray(&mymol
, 1));
572 print_rtp(ftp2fn(efRTP
, NFILE
, fnm
), "Generated by x2top",
573 atoms
, plist
, &atypes
, cgnr
);
578 dump_hybridization(debug
, atoms
, nbonds
);
580 close_symtab(&symtab
);
582 printf("\nWARNING: topologies generated by %s can not be trusted at face value.\n",
583 output_env_get_program_display_name(oenv
));
584 printf(" Please verify atomtypes and charges by comparison to other\n");
585 printf(" topologies.\n");