SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / topology / atoms.h
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) 2012,2014,2015,2016,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 #ifndef GMX_TOPOLOGY_ATOMS_H
39 #define GMX_TOPOLOGY_ATOMS_H
40
41 #include <stdio.h>
42
43 #include <vector>
44
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/real.h"
47 #include "gromacs/utility/unique_cptr.h"
48
49 struct t_symtab;
50
51 /* The particle type */
52 enum class ParticleType : int
53 {
54     Atom,
55     Nucleus,
56     Shell,
57     Bond,
58     VSite,
59     Count
60 };
61
62 /* The particle type names */
63 const char* enumValueToString(ParticleType enumValue);
64
65 /* Enumerated type for pdb records. The other entries are ignored
66  * when reading a pdb file
67  */
68 enum class PdbRecordType : int
69 {
70     Atom,
71     Hetatm,
72     Anisou,
73     Cryst1,
74     Compound,
75     Model,
76     EndModel,
77     Ter,
78     Header,
79     Title,
80     Remark,
81     Conect,
82     Count
83 };
84
85 const char* enumValueToString(PdbRecordType enumValue);
86
87 typedef struct t_atom
88 {
89     real           m, q;       /* Mass and charge                      */
90     real           mB, qB;     /* Mass and charge for Free Energy calc */
91     unsigned short type;       /* Atom type                            */
92     unsigned short typeB;      /* Atom type for Free Energy calc       */
93     ParticleType   ptype;      /* Particle type                        */
94     int            resind;     /* Index into resinfo (in t_atoms)      */
95     int            atomnumber; /* Atomic Number or 0                   */
96     char           elem[4];    /* Element name                         */
97 } t_atom;
98
99 typedef struct t_resinfo
100 {
101     char**        name;     /* Pointer to the residue name          */
102     int           nr;       /* Residue number                       */
103     unsigned char ic;       /* Code for insertion of residues       */
104     int           chainnum; /* Iincremented at TER or new chain id  */
105     char          chainid;  /* Chain identifier written/read to pdb */
106     char**        rtp;      /* rtp building block name (optional)   */
107 } t_resinfo;
108
109 typedef struct t_pdbinfo
110 {
111     PdbRecordType type;         /* PDB record name                      */
112     int           atomnr;       /* PDB atom number                      */
113     char          altloc;       /* Alternate location indicator         */
114     char          atomnm[6];    /* True atom name including leading spaces */
115     real          occup;        /* Occupancy                            */
116     real          bfac;         /* B-factor                             */
117     gmx_bool      bAnisotropic; /* (an)isotropic switch                 */
118     int           uij[6];       /* Anisotropic B-factor                 */
119 } t_pdbinfo;
120
121 //! Contains indices into group names for different groups.
122 using AtomGroupIndices = std::vector<int>;
123
124 typedef struct t_atoms
125 {
126     int     nr;         /* Nr of atoms                          */
127     t_atom* atom;       /* Array of atoms (dim: nr)             */
128                         /* The following entries will not       */
129                         /* always be used (nres==0)             */
130     char*** atomname;   /* Array of pointers to atom name       */
131                         /* use: (*(atomname[i]))                */
132     char*** atomtype;   /* Array of pointers to atom types      */
133                         /* use: (*(atomtype[i]))                */
134     char*** atomtypeB;  /* Array of pointers to B atom types    */
135                         /* use: (*(atomtypeB[i]))               */
136     int        nres;    /* The number of resinfo entries        */
137     t_resinfo* resinfo; /* Array of residue names and numbers   */
138     t_pdbinfo* pdbinfo; /* PDB Information, such as aniso. Bfac */
139
140     /* Flags that tell if properties are set for all nr atoms.
141      * For B-state parameters, both haveBState and the mass/charge/type
142      * flag should be TRUE.
143      */
144     gmx_bool haveMass;    /* Mass available                       */
145     gmx_bool haveCharge;  /* Charge available                     */
146     gmx_bool haveType;    /* Atom type available                  */
147     gmx_bool haveBState;  /* B-state parameters available         */
148     gmx_bool havePdbInfo; /* pdbinfo available                    */
149 } t_atoms;
150
151 typedef struct t_atomtypes
152 {
153     int  nr;         /* number of atomtypes                          */
154     int* atomnumber; /* Atomic number, used for QM/MM                */
155 } t_atomtypes;
156
157 #define PERTURBED(a) (((a).mB != (a).m) || ((a).qB != (a).q) || ((a).typeB != (a).type))
158
159 void init_atom(t_atoms* at);
160 void init_atomtypes(t_atomtypes* at);
161 void done_atom(t_atoms* at);
162 void done_and_delete_atoms(t_atoms* atoms);
163 void done_atomtypes(t_atomtypes* at);
164
165 void init_t_atoms(t_atoms* atoms, int natoms, gmx_bool bPdbinfo);
166 /* allocate memory for the arrays, set nr to natoms and nres to 0
167  * set pdbinfo to NULL or allocate memory for it */
168
169 void gmx_pdbinfo_init_default(t_pdbinfo* pdbinfo);
170
171 t_atoms* copy_t_atoms(const t_atoms* src);
172 /* copy an atoms struct from src to a new one */
173
174 void add_t_atoms(t_atoms* atoms, int natom_extra, int nres_extra);
175 /* allocate extra space for more atoms and or residues */
176
177 void t_atoms_set_resinfo(t_atoms*         atoms,
178                          int              atom_ind,
179                          struct t_symtab* symtab,
180                          const char*      resname,
181                          int              resnr,
182                          unsigned char    ic,
183                          int              chainnum,
184                          char             chainid);
185 /* Set the residue name, number, insertion code and chain identifier
186  * of atom index atom_ind.
187  */
188
189 void pr_atoms(FILE* fp, int indent, const char* title, const t_atoms* atoms, gmx_bool bShownumbers);
190 void pr_atomtypes(FILE* fp, int indent, const char* title, const t_atomtypes* atomtypes, gmx_bool bShowNumbers);
191
192 /*! \brief Compare information in the t_atoms data structure.
193  *
194  * \param[in] fp Pointer to file to write to.
195  * \param[in] a1 Pointer to first data structure to compare.
196  * \param[in] a2 Pointer to second data structure or nullptr.
197  * \param[in] relativeTolerance Relative floating point comparison tolerance.
198  * \param[in] absoluteTolerance Absolute floating point comparison tolerance.
199  */
200 void compareAtoms(FILE* fp, const t_atoms* a1, const t_atoms* a2, real relativeTolerance, real absoluteTolerance);
201
202 /*! \brief Set mass for each atom using the atom and residue names using a database
203  *
204  * If atoms->haveMass = TRUE does nothing.
205  * If printMissingMasss = TRUE, prints details for first 10 missing masses
206  * to stderr.
207  */
208 void atomsSetMassesBasedOnNames(t_atoms* atoms, gmx_bool printMissingMasses);
209
210 //! Deleter for t_atoms, needed until it has a proper destructor.
211 using AtomsDataPtr = gmx::unique_cptr<t_atoms, done_and_delete_atoms>;
212
213
214 #endif