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