Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / fileio / pdbio.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) 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37
38 #ifndef GMX_FILEIO_PDBIO_H
39 #define GMX_FILEIO_PDBIO_H
40
41 #include <stdio.h>
42
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/topology/atoms.h"
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/real.h"
47
48 class AtomProperties;
49 struct t_atoms;
50 struct t_symtab;
51 struct t_topology;
52
53 typedef struct gmx_conect_t* gmx_conect;
54
55 /* Write a PDB line with an ATOM or HETATM record directly to a file pointer.
56  *
57  * Returns the number of characters printed.
58  */
59 int gmx_fprintf_pdb_atomline(FILE*           fp,
60                              enum PDB_record record,
61                              int             atom_seq_number,
62                              const char*     atom_name,
63                              char            alternate_location,
64                              const char*     res_name,
65                              char            chain_id,
66                              int             res_seq_number,
67                              char            res_insertion_code,
68                              real            x,
69                              real            y,
70                              real            z,
71                              real            occupancy,
72                              real            b_factor,
73                              const char*     element);
74
75 /* Enumerated value for indexing an uij entry (anisotropic temperature factors) */
76 enum
77 {
78     U11,
79     U22,
80     U33,
81     U12,
82     U13,
83     U23
84 };
85
86 void pdb_use_ter(gmx_bool bSet);
87 /* set read_pdbatoms to read upto 'TER' or 'ENDMDL' (default, bSet=FALSE).
88    This function is fundamentally broken as far as thread-safety is concerned.*/
89
90 void gmx_write_pdb_box(FILE* out, int ePBC, const matrix box);
91 /* write the box in the CRYST1 record,
92  * with ePBC=-1 the pbc is guessed from the box
93  * This function is fundamentally broken as far as thread-safety is concerned.
94  */
95
96 void write_pdbfile_indexed(FILE*          out,
97                            const char*    title,
98                            const t_atoms* atoms,
99                            const rvec     x[],
100                            int            ePBC,
101                            const matrix   box,
102                            char           chain,
103                            int            model_nr,
104                            int            nindex,
105                            const int      index[],
106                            gmx_conect     conect,
107                            bool           usePqrFormat);
108 /* REALLY low level */
109
110 void write_pdbfile(FILE*          out,
111                    const char*    title,
112                    const t_atoms* atoms,
113                    const rvec     x[],
114                    int            ePBC,
115                    const matrix   box,
116                    char           chain,
117                    int            model_nr,
118                    gmx_conect     conect);
119 /* Low level pdb file writing routine.
120  *
121  *          ONLY FOR SPECIAL PURPOSES,
122  *
123  *       USE write_sto_conf WHEN YOU CAN.
124  *
125  * override chain-identifiers with chain when chain>0
126  * write ENDMDL when bEndmodel is TRUE.
127  *
128  * If the gmx_conect structure is not NULL its content is dumped as CONECT records
129  * which may be useful for visualization purposes.
130  */
131
132 void get_pdb_atomnumber(const t_atoms* atoms, AtomProperties* aps);
133 /* Routine to extract atomic numbers from the atom names */
134
135 int read_pdbfile(FILE*            in,
136                  char*            title,
137                  int*             model_nr,
138                  struct t_atoms*  atoms,
139                  struct t_symtab* symtab,
140                  rvec             x[],
141                  int*             ePBC,
142                  matrix           box,
143                  gmx_bool         bChange,
144                  gmx_conect       conect);
145 /* Function returns number of atoms found.
146  * ePBC and gmx_conect structure may be NULL.
147  */
148
149 void gmx_pdb_read_conf(const char* infile, t_symtab* symtab, char** name, t_atoms* atoms, rvec x[], int* ePBC, matrix box);
150 /* Read a pdb file and extract ATOM and HETATM fields.
151  * Read a box from the CRYST1 line, return 0 box when no CRYST1 is found.
152  * ePBC may be NULL.
153  *
154  * If name is not nullptr, gmx_strdup the title string into it. */
155
156 void get_pdb_coordnum(FILE* in, int* natoms);
157 /* Read a pdb file and count the ATOM and HETATM fields. */
158
159 gmx_bool is_hydrogen(const char* nm);
160 /* Return whether atom nm is a hydrogen */
161
162 gmx_bool is_dummymass(const char* nm);
163 /* Return whether atom nm is a dummy mass */
164
165 /* Routines to handle CONECT records if they have been read in */
166 void gmx_conect_dump(FILE* fp, gmx_conect conect);
167
168 gmx_bool gmx_conect_exist(gmx_conect conect, int ai, int aj);
169 /* Return TRUE if there is a conection between the atoms */
170
171 void gmx_conect_add(gmx_conect conect, int ai, int aj);
172 /* Add a connection between ai and aj (numbered from 0 to natom-1) */
173
174 gmx_conect gmx_conect_generate(const t_topology* top);
175 /* Generate a conect structure from a topology */
176
177 gmx_conect gmx_conect_init();
178 /* Initiate data structure */
179
180 void gmx_conect_done(gmx_conect gc);
181 /* Free memory */
182
183 #endif /* GMX_FILEIO_PDBIO_H */