Fixing copyright issues and code contributors
[alexxy/gromacs.git] / src / gmxlib / molfile_plugin.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * This file is part of Gromacs        Copyright (c) 1991-2008
5  * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
6  * Copyright (c) 2012,2013, by the GROMACS development team, led by
7  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
8  * others, as listed in the AUTHORS file in the top-level source
9  * 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 /***************************************************************************
39  *cr
40  *cr            (C) Copyright 1995-2006 The Board of Trustees of the
41  *cr                        University of Illinois
42  *cr                         All Rights Reserved
43  *cr
44
45 Developed by:           Theoretical and Computational Biophysics Group
46                         University of Illinois at Urbana-Champaign
47                         http://www.ks.uiuc.edu/
48
49 Permission is hereby granted, free of charge, to any person obtaining a copy of
50 this software and associated documentation files (the Software), to deal with
51 the Software without restriction, including without limitation the rights to
52 use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies
53 of the Software, and to permit persons to whom the Software is furnished to
54 do so, subject to the following conditions:
55
56 Redistributions of source code must retain the above copyright notice,
57 this list of conditions and the following disclaimers.
58
59 Redistributions in binary form must reproduce the above copyright notice,
60 this list of conditions and the following disclaimers in the documentation
61 and/or other materials provided with the distribution.
62
63 Neither the names of Theoretical and Computational Biophysics Group,
64 University of Illinois at Urbana-Champaign, nor the names of its contributors
65 may be used to endorse or promote products derived from this Software without
66 specific prior written permission.
67
68 THE SOFTWARE IS PROVIDED AS IS, WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
69 IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
70 FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.  IN NO EVENT SHALL
71 THE CONTRIBUTORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR
72 OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE,
73 ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR
74 OTHER DEALINGS WITH THE SOFTWARE.
75  ***************************************************************************/
76
77 /***************************************************************************
78  * RCS INFORMATION:
79  *
80  *      $RCSfile: molfile_plugin.h,v $
81  *      $Author: saam $       $Locker:  $             $State: Exp $
82  *      $Revision: 1.91 $       $Date: 2009/07/28 21:54:30 $
83  *
84  ***************************************************************************/
85
86 /** @file 
87  * API for C extensions to define a way to load structure, coordinate,
88  * trajectory, and volumetric data files  
89  */ 
90
91 #ifndef MOL_FILE_PLUGIN_H
92 #define MOL_FILE_PLUGIN_H
93
94 #include "vmdplugin.h"
95
96 /**
97  * Define a common plugin type to be used when registering the plugin.
98  */
99 #define MOLFILE_PLUGIN_TYPE "mol file reader"
100
101 /**
102  * File converter plugins use the same API  but register under a different
103  * type so that regular file readers can have priority.
104  */
105 #define MOLFILE_CONVERTER_PLUGIN_TYPE "mol file converter"
106
107 /* File plugin symbolic constants for better code readability */
108 #define MOLFILE_SUCCESS           0   /**< succeeded in reading file      */
109 #define MOLFILE_EOF              -1   /**< end of file                    */
110 #define MOLFILE_ERROR            -1   /**< error reading/opening a file   */
111 #define MOLFILE_NOSTRUCTUREDATA  -2   /**< no structure data in this file */
112
113 #define MOLFILE_NUMATOMS_UNKNOWN -1   /**< unknown number of atoms       */
114 #define MOLFILE_NUMATOMS_NONE     0   /**< no atoms in this file type    */
115
116 /**
117  * Maximum string size macro
118  */
119 #define MOLFILE_BUFSIZ           81   /**< maximum chars in string data  */
120 #define MOLFILE_BIGBUFSIZ        4096 /**< maximum chars in long strings */
121
122 #define MOLFILE_MAXWAVEPERTS     25   /**< maximum number of wavefunctions
123                                        *   per timestep */
124
125 #define MOLFILE_QM_STATUS_UNKNOWN  -1
126 #define MOLFILE_QM_OPT_CONVERGED    0
127 #define MOLFILE_QM_SCF_NOT_CONV     1
128 #define MOLFILE_QM_OPT_NOT_CONV     2
129 #define MOLFILE_QM_FILE_TRUNCATED   3
130
131
132 /**
133  * File level comments, origin information, and annotations.
134  */
135 typedef struct {
136   char database[81];   /**< database of origin, if any        */
137   char accession[81];  /**< database accession code, if any   */
138   char date[81];       /**< date/time stamp for this data     */
139   char title[81];      /**< brief title for this data         */
140   int remarklen;       /**< length of remarks string          */
141   char *remarks;       /**< free-form remarks about data      */
142 } molfile_metadata_t;
143
144
145 /* 
146  * Struct for specifying atoms in a molecular structure.  The first 
147  * six components are required, the rest are optional and their presence is 
148  * indicating by setting the corresponding bit in optsflag.  When omitted,
149  * the application (for read_structure) or plugin (for write_structure) 
150  * must be able to supply default values if the missing parameters are 
151  * part of its internal data structure.
152  * Note that it is not possible to specify coordinates with this structure.
153  * This is intentional; all coordinate I/O is done with the read_timestep and 
154  * write_timestep functions. 
155  */
156
157 /**
158  * Per-atom attributes and information.
159  */
160 typedef struct {
161   /* these fields absolutely must be set or initialized to empty */
162   char name[16];      /**< required atom name string             */
163   char type[16];      /**< required atom type string             */
164   char resname[8];    /**< required residue name string          */
165   int resid;          /**< required integer residue ID           */
166   char segid[8];      /**< required segment name string, or ""   */
167   char chain[2];      /**< required chain name, or ""            */
168
169   /* rest are optional; use optflags to specify what's present   */
170   char altloc[2];     /**< optional PDB alternate location code  */
171   char insertion[2];  /**< optional PDB insertion code           */
172   float occupancy;    /**< optional occupancy value              */
173   float bfactor;      /**< optional B-factor value               */
174   float mass;         /**< optional mass value                   */
175   float charge;       /**< optional charge value                 */
176   float radius;       /**< optional radius value                 */
177   int atomicnumber;   /**< optional element atomic number        */
178 } molfile_atom_t;
179
180 /*@{*/
181 /** Plugin optional data field availability flag */
182 #define MOLFILE_NOOPTIONS     0x0000 /**< no optional data                 */
183 #define MOLFILE_INSERTION     0x0001 /**< insertion codes provided         */
184 #define MOLFILE_OCCUPANCY     0x0002 /**< occupancy data provided          */
185 #define MOLFILE_BFACTOR       0x0004 /**< B-factor data provided           */
186 #define MOLFILE_MASS          0x0008 /**< Atomic mass provided             */
187 #define MOLFILE_CHARGE        0x0010 /**< Atomic charge provided           */
188 #define MOLFILE_RADIUS        0x0020 /**< Atomic VDW radius provided       */
189 #define MOLFILE_ALTLOC        0x0040 /**< Multiple conformations present   */
190 #define MOLFILE_ATOMICNUMBER  0x0080 /**< Atomic element number provided   */
191 #define MOLFILE_BONDSSPECIAL  0x0100 /**< Only non-standard bonds provided */
192 #define MOLFILE_BADOPTIONS    0xFFFFFFFF /**< Detect badly behaved plugins */
193                               
194 /*@}*/
195
196 /*@{*/
197 /** Plugin optional data field availability flag */
198 #define MOLFILE_QMTS_NOOPTIONS     0x0000 /**< no optional data               */
199 #define MOLFILE_QMTS_GRADIENT      0x0001 /**< energy gradients provided      */
200 #define MOLFILE_QMTS_SCFITER       0x0002
201 /*@}*/
202
203 #if vmdplugin_ABIVERSION > 10
204 typedef struct molfile_timestep_metadata {
205   unsigned int count;                  /**< total # timesteps; -1 if unknown */
206   unsigned int avg_bytes_per_timestep; /** bytes per timestep                */
207   int has_velocities;                  /**< if timesteps have velocities     */
208 } molfile_timestep_metadata_t;
209 #endif
210
211 #if vmdplugin_ABIVERSION > 11
212 typedef struct molfile_qm_timestep_metadata {
213   unsigned int count;                  /**< total # timesteps; -1 if unknown */
214   unsigned int avg_bytes_per_timestep; /** bytes per timestep                */
215   int has_gradient;                    /**< if timestep contains gradient    */
216   int num_scfiter;                     /**< # scf iterations for this ts     */
217   int num_orbitals_per_wavef[MOLFILE_MAXWAVEPERTS]; /**< # orbitals for each wavefunction */
218   int has_orben_per_wavef[MOLFILE_MAXWAVEPERTS]; /**< orbital energy flags */
219   int has_occup_per_wavef[MOLFILE_MAXWAVEPERTS]; /**< orbital occupancy flags */
220   int num_wavef ;                      /**< # wavefunctions in this ts     */
221   int wavef_size;                      /**< size of one wavefunction 
222                                         *   (# of gaussian basis fctns)    */
223   int num_charge_sets;            /**< # of charge values per atom */
224 } molfile_qm_timestep_metadata_t;
225 #endif
226
227 /*
228  * Per-timestep atom coordinates and periodic cell information
229  */ 
230 typedef struct {
231   float *coords;        /**< coordinates of all atoms, arranged xyzxyzxyz   */
232 #if vmdplugin_ABIVERSION > 10
233   float *velocities;    /**< space for velocities of all atoms; same layout */
234                         /**< NULL unless has_velocities is set              */
235 #endif
236
237   /*@{*/   
238   /**
239    * Unit cell specification of the form A, B, C, alpha, beta, gamma.
240    * notes: A, B, C are side lengths of the unit cell
241    * alpha = angle between b and c
242    *  beta = angle between a and c
243    * gamma = angle between a and b
244    */ 
245   float A, B, C, alpha, beta, gamma; 
246   /*@}*/   
247
248 #if vmdplugin_ABIVERSION > 10
249   double physical_time; /**< physical time point associated with this frame */
250 #endif
251 } molfile_timestep_t;
252
253
254 /**
255  * Metadata for volumetric datasets, read initially and used for subsequent
256  * memory allocations and file loading.  
257  */
258 typedef struct {
259   char dataname[256];   /**< name of volumetric data set                    */
260   float origin[3];      /**< origin: origin of volume (x=0, y=0, z=0 corner */
261
262   /*
263    * x/y/z axis:
264    * These the three cell sides, providing both direction and length
265    * (not unit vectors) for the x, y, and z axes.  In the simplest
266    * case, these would be <size,0,0> <0,size,0> and <0,0,size) for 
267    * an orthogonal cubic volume set.  For other cell shapes these
268    * axes can be oriented non-orthogonally, and the parallelpiped
269    * may have different side lengths, not just a cube/rhombus.
270    */
271   float xaxis[3];       /**< direction (and length) for X axis              */ 
272   float yaxis[3];       /**< direction (and length) for Y axis              */
273   float zaxis[3];       /**< direction (and length) for Z axis              */
274
275   /*
276    * x/y/z size: 
277    * Number of grid cells along each axis.  This is _not_ the
278    * physical size of the box, this is the number of voxels in each
279    * direction, independent of the shape of the volume set. 
280    */
281   int xsize;            /**< number of grid cells along the X axis          */
282   int ysize;            /**< number of grid cells along the Y axis          */
283   int zsize;            /**< number of grid cells along the Z axis          */
284
285   int has_color;        /**< flag indicating presence of voxel color data   */
286 } molfile_volumetric_t;
287
288
289 #if vmdplugin_ABIVERSION > 9
290
291 /**
292  * Sizes of various QM-related data arrays which must be allocated by
293  * the caller (VMD) so that the plugin can fill in the arrays with data.
294  */
295 typedef struct {
296   /* hessian data */
297   int nimag;                    /**< number of imaginary modes */
298   int nintcoords;               /**< number internal coordinates */
299   int ncart;                    /**< number cartesian coordinates */
300
301   /* orbital/basisset data */
302   int num_basis_funcs;          /**< number of uncontracted basis functions in basis array */
303   int num_basis_atoms;          /**< number of atoms in basis set */
304   int num_shells;               /**< total number of atomic shells */
305   int wavef_size;               /**< size of the wavefunction
306                                  *   i.e. size of secular eq. or
307                                  *   # of cartesian contracted
308                                  *   gaussian basis functions */
309
310   /* everything else */
311   int have_sysinfo;
312   int have_carthessian;         /**< hessian in cartesian coords available  */
313   int have_inthessian;          /**< hessian in internal coords available  */
314   int have_normalmodes;         /**< normal modes available  */
315 } molfile_qm_metadata_t;
316
317
318 /**
319  * struct holding the data of hessian/normal mode runs
320  * needed to calculate bond/angle/dihedral force constants
321  * XXX: do we really need doubles here??
322  */
323 typedef struct {
324   double *carthessian;  /**< hessian matrix in cartesian coordinates (ncart)*(ncart)
325                          *   as a single array of doubles (row(1), ...,row(natoms)) */
326   int    *imag_modes;   /**< list(nimag) of imaginary modes */
327   double *inthessian;   /**< hessian matrix in internal coordinates
328                          *   (nintcoords*nintcoords) as a single array of
329                          *   doubles (row(1), ...,row(nintcoords)) */
330   float *wavenumbers;   /**< array(ncart) of wavenumbers of normal modes */
331   float *intensities;   /**< array(ncart) of intensities of normal modes */
332   float *normalmodes;   /**< matrix(ncart*ncart) of normal modes  */
333 } molfile_qm_hessian_t;
334
335
336 /**
337  * struct holding the data for wavefunction/orbitals
338  * needed to generate the volumetric orbital data
339  */
340 typedef struct {
341   int *num_shells_per_atom; /**< number of shells per atom */
342   int *num_prim_per_shell;  /**< number of shell primitives shell */
343
344   float *basis;             /**< contraction coeffients and exponents for
345                              *   the basis functions in the form
346                              *   { exp(1), c-coeff(1), exp(2), c-coeff(2), ....};
347                              *   size=2*num_basis_funcs */
348   int *atomic_number;       /**< atomic numbers (chem. element) of atoms in basis set */
349   int *angular_momentum;    /**< 3 ints per wave function coefficient do describe the 
350                              *   cartesian components of the angular momentum.
351                              *   E.g. S={0 0 0}, Px={1 0 0}, Dxy={1 1 0}, or Fyyz={0 2 1}. 
352                              */
353   int *shell_symmetry;      /**< symmetry type per shell in basis */
354 } molfile_qm_basis_t;
355
356
357 /**
358  * QM run info. Parameters that stay unchanged during a single file.
359  */ 
360 typedef struct {
361   int nproc;             /**< number of processors used. XXX:? */
362   int memory;            /**< amount of memory used in Mbyte. XXX:? */ 
363   int runtype;           /**< flag indicating the calculation method. */
364   int scftype;           /**< SCF type: RHF, UHF, ROHF, GVB or MCSCF wfn. */
365   int status;            /**< indicates wether SCF and geometry optimization
366                           *   have converged properly. */
367   int num_electrons;     /**< number of electrons.    XXX: can be fractional in some DFT codes */
368   int totalcharge;       /**< total charge of system. XXX: can be fractional in some DFT codes */
369   int num_occupied_A;    /**< number of occupied alpha orbitals */
370   int num_occupied_B;    /**< number of occupied beta orbitals */
371
372   double *nuc_charge;    /**< array(natom) containing the nuclear charge of atom i */
373
374   char basis_string[MOLFILE_BUFSIZ];    /**< basis name as "nice" string. */
375   char runtitle[MOLFILE_BIGBUFSIZ];     /**< title of run.                */
376   char geometry[MOLFILE_BUFSIZ];        /**< type of provided geometry,   XXX: remove?
377                                          * e.g. UNIQUE, ZMT, CART, ...    */
378   char version_string[MOLFILE_BUFSIZ];  /**< QM code version information. */
379 } molfile_qm_sysinfo_t;
380
381
382 typedef struct {
383   int   type;               /**< CANONICAL, LOCALIZED, OTHER */
384   int   spin;               /**< 1 for alpha, -1 for beta */
385   int   excitation;         /**< 0 for ground state, 1,2,3,... for excited states */
386   int   multiplicity;       /**< spin multiplicity of the state, zero if unknown */
387   char info[MOLFILE_BUFSIZ]; /**< string for additional type info */
388
389   double energy;            /**< energy of the electronic state.
390                              *   i.e. HF-SCF energy, CI state energy,
391                              *   MCSCF energy, etc. */
392
393   float *wave_coeffs;       /**< expansion coefficients for wavefunction in the
394                              *   form {orbital1(c1),orbital1(c2),.....,orbitalM(cN)} */
395   float *orbital_energies;  /**< list of orbital energies for wavefunction */
396   float *occupancies;       /**< orbital occupancies */
397   int   *orbital_ids;       /**< orbital ID numbers; If NULL then VMD will
398                              *   assume 1,2,3,...num_orbs.     */
399 } molfile_qm_wavefunction_t;
400
401
402 /**
403  * QM per trajectory timestep info
404  */
405 typedef struct {
406   molfile_qm_wavefunction_t *wave; /**< array of wavefunction objects */
407   float  *gradient;         /**< force on each atom (=gradient of energy) */
408
409   double *scfenergies;      /**< scfenergy per trajectory point. */
410   double *charges;          /**< per-atom charges */
411   int    *charge_types;     /**< type of each charge set */
412 } molfile_qm_timestep_t;
413
414
415 /**
416  * QM wavefunctions, and related information 
417  */
418 typedef struct {
419   molfile_qm_hessian_t hess;            /* hessian info */
420   molfile_qm_basis_t   basis;           /* basis set info */
421   molfile_qm_sysinfo_t run;             /* system info  */
422 } molfile_qm_t;
423
424
425 /**
426  *  Enumeration of all of the wavefunction types that can be read
427  *  from QM file reader plugins.
428  *
429  * CANON    = canonical (i.e diagonalized) wavefunction
430  * GEMINAL  = GVB-ROHF geminal pairs
431  * MCSCFNAT = Multi-Configuration SCF natural orbitals
432  * MCSCFOPT = Multi-Configuration SCF optimized orbitals
433  * CINATUR  = Configuration-Interaction natural orbitals
434  * BOYS     = Boys localization
435  * RUEDEN   = Ruedenberg localization
436  * PIPEK    = Pipek-Mezey population localization
437  *
438  * NBO related localizations:
439  * --------------------------
440  * NAO      = Natural Atomic Orbitals
441  * PNAO     = pre-orthogonal NAOs
442  * NBO      = Natural Bond Orbitals
443  * PNBO     = pre-orthogonal NBOs
444  * NHO      = Natural Hybrid Orbitals
445  * PNHO     = pre-orthogonal NHOs
446  * NLMO     = Natural Localized Molecular Orbitals
447  * PNLMO    = pre-orthogonal NLMOs
448  *
449  * UNKNOWN  = Use this for any type not listed here
450  *            You can use the string field for description
451  */
452 enum molfile_qm_wavefunc_type {
453   MOLFILE_WAVE_CANON,    MOLFILE_WAVE_GEMINAL,
454   MOLFILE_WAVE_MCSCFNAT, MOLFILE_WAVE_MCSCFOPT,
455   MOLFILE_WAVE_CINATUR,
456   MOLFILE_WAVE_PIPEK,  MOLFILE_WAVE_BOYS, MOLFILE_WAVE_RUEDEN,
457   MOLFILE_WAVE_NAO,    MOLFILE_WAVE_PNAO, MOLFILE_WAVE_NHO, 
458   MOLFILE_WAVE_PNHO,   MOLFILE_WAVE_NBO,  MOLFILE_WAVE_PNBO, 
459   MOLFILE_WAVE_PNLMO,  MOLFILE_WAVE_NLMO, MOLFILE_WAVE_MOAO, 
460   MOLFILE_WAVE_NATO,   MOLFILE_WAVE_UNKNOWN
461 };
462
463 enum molfile_qm_charge_type {
464   MOLFILE_QMCHARGE_MULLIKEN, MOLFILE_QMCHARGE_LOWDIN,
465   MOLFILE_QMCHARGE_ESP, MOLFILE_QMCHARGE_NPA,
466   MOLFILE_QMCHARGE_UNKNOWN
467 };
468 #endif
469
470
471 /**
472  *  Enumeration of all of the supported graphics objects that can be read
473  *  from graphics file reader plugins.
474  */
475 enum molfile_graphics_type {
476   MOLFILE_POINT,  MOLFILE_TRIANGLE, MOLFILE_TRINORM, MOLFILE_NORMS, 
477   MOLFILE_LINE,   MOLFILE_CYLINDER, MOLFILE_CAPCYL,  MOLFILE_CONE,    
478   MOLFILE_SPHERE, MOLFILE_TEXT,     MOLFILE_COLOR,   MOLFILE_TRICOLOR
479 };
480
481 /**
482  *  Individual graphics object/element data
483  */ 
484 typedef struct {
485   int type;             /* One of molfile_graphics_type */
486   int style;            /* A general style parameter    */
487   float size;           /* A general size parameter     */
488   float data[9];        /* All data for the element     */
489 } molfile_graphics_t;
490
491
492 /*
493  * Types for raw graphics elements stored in files.  Data for each type
494  * should be stored by the plugin as follows:
495
496 type        data                                     style       size
497 ----        ----                                     -----       ----
498 point       x, y, z                                              pixel size
499 triangle    x1,y1,z1,x2,y2,z2,x3,y3,z3                 
500 trinorm     x1,y1,z1,x2,y2,z2,x3,y3,z3                 
501             the next array element must be NORMS
502 tricolor    x1,y1,z1,x2,y2,z2,x3,y3,z3                 
503             the next array elements must be NORMS
504             the following element must be COLOR, with three RGB triples
505 norms       x1,y1,z1,x2,y2,z2,x3,y3,z3                 
506 line        x1,y1,z1,x2,y2,z2                        0=solid     pixel width
507                                                      1=stippled
508 cylinder    x1,y1,z1,x2,y2,z2                        resolution  radius
509 capcyl      x1,y1,z1,x2,y2,z2                        resolution  radius
510 sphere      x1,y1,z1                                 resolution  radius
511 text        x, y, z, up to 24 bytes of text                      pixel size
512 color       r, g, b
513 */
514
515
516 /**
517  * Main file reader API.  Any function in this struct may be NULL
518  * if not implemented by the plugin; the application checks this to determine
519  * what functionality is present in the plugin. 
520  */ 
521 typedef struct {
522   /**
523    * Required header 
524    */
525   vmdplugin_HEAD;
526
527   /**
528    * Filename extension for this file type.  May be NULL if no filename 
529    * extension exists and/or is known.  For file types that match several
530    * common extensions, list them in a comma separated list such as:
531    *  "pdb,ent,foo,bar,baz,ban"
532    * The comma separated list will be expanded when filename extension matching
533    * is performed.  If multiple plugins solicit the same filename extensions,
534    * the one that lists the extension earliest in its list is selected. In the 
535    * case of a "tie", the first one tried/checked "wins".
536    */
537   const char *filename_extension;
538
539   /**
540    * Try to open the file for reading.  Return an opaque handle, or NULL on
541    * failure. Set the number of atoms; if the number of atoms cannot be 
542    * determined, set natoms to MOLFILE_NUMATOMS_UNKNOWN. 
543    * Filetype should be the name under which this plugin was registered;
544    * this is provided so that plugins can provide the same function pointer
545    * to handle multiple file types.
546    */
547   void *(* open_file_read)(const char *filepath, const char *filetype, 
548       int *natoms);
549   
550   /**
551    * Read molecular structure from the given file handle.  atoms is allocated
552    * by the caller and points to space for natoms.
553    * On success, place atom information in the passed-in pointer.  
554    * optflags specifies which optional fields in the atoms will be set by
555    * the plugin.
556    */
557   int (*read_structure)(void *, int *optflags, molfile_atom_t *atoms);
558
559   /**
560    * Read bond information for the molecule.  On success the arrays from
561    * and to should point to the (one-based) indices of bonded atoms.
562    * Each unique bond should be specified only once, so file formats that list
563    * bonds twice will need post-processing before the results are returned to
564    * the caller.
565    * If the plugin provides bond information, but the file loaded doesn't 
566    * actually contain any bond info, the nbonds parameter should be
567    * set to 0 and from/to should be set to NULL to indicate that no bond
568    * information was actually present, and automatic bond search should be
569    * performed.  
570    *
571    * If the plugin provides bond order information, the bondorder array
572    * will contain the bond order for each from/to pair.  If not, the bondorder
573    * pointer should be set to NULL, in which case the caller will provide a 
574    * default bond order value of 1.0.
575    *
576    * If the plugin provides bond type information, the bondtype array
577    * will contain the bond type index for each from/to pair. These numbers
578    * are consecutive integers starting from 0.
579    * the bondtypenames list, contains the corresponding names, if available,
580    * as a NULL string terminated list. nbondtypes is provided for convenience
581    * and consistency checking.
582    *
583    * These arrays must be freed by the plugin in the close_file_read function.
584    * This function can be called only after read_structure().  
585    * Return MOLFILE_SUCCESS if no errors occur. 
586    */
587 #if vmdplugin_ABIVERSION > 14
588   int (*read_bonds)(void *, int *nbonds, int **from, int **to, float **bondorder, 
589                     int **bondtype, int *nbondtypes, char ***bondtypename);
590 #else
591   int (*read_bonds)(void *, int *nbonds, int **from, int **to, float **bondorder);
592 #endif
593
594   /**
595    * XXX this function will be augmented and possibly superceded by a 
596    *     new QM-capable version named read_timestep(), when finished.
597    *
598    * Read the next timestep from the file.  Return MOLFILE_SUCCESS, or 
599    * MOLFILE_EOF on EOF.  If the molfile_timestep_t argument is NULL, then 
600    * the frame should be skipped.  Otherwise, the application must prepare 
601    * molfile_timestep_t by allocating space in coords for the corresponding 
602    * number of coordinates.  
603    * The natoms parameter exists because some coordinate file formats 
604    * (like CRD) cannot determine for themselves how many atoms are in a 
605    * timestep; the app must therefore obtain this information elsewhere
606    * and provide it to the plugin.
607    */
608   int (* read_next_timestep)(void *, int natoms, molfile_timestep_t *);
609
610   /**
611    * Close the file and release all data.  The handle cannot be reused.
612    */
613   void (* close_file_read)(void *);
614    
615   /**
616    * Open a coordinate file for writing using the given header information.
617    * Return an opaque handle, or NULL on failure.  The application must
618    * specify the number of atoms to be written. 
619    * filetype should be the name under which this plugin was registered.
620    */
621   void *(* open_file_write)(const char *filepath, const char *filetype, 
622       int natoms);
623   
624   /**
625    * Write structure information.  Return success.
626    */
627   int (* write_structure)(void *, int optflags, const molfile_atom_t *atoms);
628
629   /**
630    * Write a timestep to the coordinate file.  Return MOLFILE_SUCCESS if no
631    * errors occur.  If the file contains structure information in each 
632    * timestep (like a multi-entry PDB), it will have to cache the information 
633    * from the initial calls from write_structure.
634    */
635   int (* write_timestep)(void *, const molfile_timestep_t *);
636   
637   /**
638    * Close the file and release all data.  The handle cannot be reused.
639    */
640   void (* close_file_write)(void *);
641
642   /**
643    * Retrieve metadata pertaining to volumetric datasets in this file.
644    * Set nsets to the number of volumetric data sets, and set *metadata
645    * to point to an array of molfile_volumetric_t.  The array is owned by
646    * the plugin and should be freed by close_file_read().  The application
647    * may call this function any number of times.
648    */
649   int (* read_volumetric_metadata)(void *, int *nsets, 
650         molfile_volumetric_t **metadata);
651
652   /** 
653    * Read the specified volumetric data set into the space pointed to by 
654    * datablock.  The set is specified with a zero-based index.  The space 
655    * allocated for the datablock must be equal to
656    * xsize * ysize * zsize.  No space will be allocated for colorblock 
657    * unless has_color is nonzero; in that case, colorblock should be
658    * filled in with three RGB floats per datapoint.
659    */
660   int (* read_volumetric_data)(void *, int set, float *datablock, 
661         float *colorblock);
662
663   /**
664    * Read raw graphics data stored in this file.   Return the number of data
665    * elements and the data itself as an array of molfile_graphics_t in the 
666    * pointer provided by the application.  The plugin is responsible for 
667    * freeing the data when the file is closed.
668    */
669   int (* read_rawgraphics)(void *, int *nelem, const molfile_graphics_t **data);
670
671   /**
672    * Read molecule metadata such as what database (if any) this file/data
673    * came from, what the accession code for the database is, textual remarks
674    * and other notes pertaining to the contained structure/trajectory/volume
675    * and anything else that's informative at the whole file level.
676    */ 
677   int (* read_molecule_metadata)(void *, molfile_metadata_t **metadata);
678   
679   /**
680    * Write bond information for the molecule.  The arrays from
681    * and to point to the (one-based) indices of bonded atoms.
682    * Each unique bond will be specified only once by the caller. 
683    * File formats that list bonds twice will need to emit both the 
684    * from/to and to/from versions of each.
685    * This function must be called before write_structure().  
686    *
687    * Like the read_bonds() routine, the bondorder pointer is set to NULL
688    * if the caller doesn't have such information, in which case the 
689    * plugin should assume a bond order of 1.0 if the file format requires
690    * bond order information.
691    *
692    * Support for bond types follows the bondorder rules. bondtype is
693    * an integer array of the size nbonds that contains the bond type
694    * index (consecutive integers starting from 0) and bondtypenames
695    * contain the corresponding strings, in case the naming/numbering
696    * scheme is different from the index numbers.
697    * if the pointers are set to NULL, then this information is not available.
698    * bondtypenames can only be used of bondtypes is also given.
699    * Return MOLFILE_SUCCESS if no errors occur. 
700    */
701 #if vmdplugin_ABIVERSION > 14
702   int (* write_bonds)(void *, int nbonds, int *from, int *to, float *bondorder, 
703                      int *bondtype, int nbondtypes, char **bondtypename);
704 #else
705   int (* write_bonds)(void *, int nbonds, int *from, int *to, float *bondorder);
706 #endif
707
708 #if vmdplugin_ABIVERSION > 9
709   /**
710    * Write the specified volumetric data set into the space pointed to by 
711    * datablock.  The * allocated for the datablock must be equal to
712    * xsize * ysize * zsize.  No space will be allocated for colorblock 
713    * unless has_color is nonzero; in that case, colorblock should be
714    * filled in with three RGB floats per datapoint.
715    */
716   int (* write_volumetric_data)(void *, molfile_volumetric_t *metadata,
717                                 float *datablock, float *colorblock);
718
719 #if vmdplugin_ABIVERSION > 15
720   /** 
721    * Read in Angles, Dihedrals, Impropers, and Cross Terms and optionally types.
722    * (Cross terms pertain to the CHARMM/NAMD CMAP feature) 
723    */
724   int (* read_angles)(void *handle, int *numangles, int **angles, int **angletypes,
725                       int *numangletypes, char ***angletypenames, int *numdihedrals,
726                       int **dihedrals, int **dihedraltypes, int *numdihedraltypes,
727                       char ***dihedraltypenames, int *numimpropers, int **impropers,        
728                       int **impropertypes, int *numimpropertypes, char ***impropertypenames,
729                       int *numcterms, int **cterms, int *ctermcols, int *ctermrows);
730
731   /** 
732    * Write out Angles, Dihedrals, Impropers, and Cross Terms
733    * (Cross terms pertain to the CHARMM/NAMD CMAP feature) 
734    */
735   int (* write_angles)(void *handle, int numangles, const int *angles, const int *angletypes,
736                        int numangletypes, const char **angletypenames, int numdihedrals,
737                        const int *dihedrals, const int *dihedraltypes, int numdihedraltypes,
738                        const char **dihedraltypenames, int numimpropers, 
739                        const int *impropers, const int *impropertypes, int numimpropertypes,
740                        const char **impropertypenames, int numcterms,  const int *cterms, 
741                        int ctermcols, int ctermrows);
742 #else
743   /** 
744    * Read in Angles, Dihedrals, Impropers, and Cross Terms
745    * Forces are in Kcal/mol
746    * (Cross terms pertain to the CHARMM/NAMD CMAP feature, forces are given
747    *  as a 2-D matrix)
748    */
749   int (* read_angles)(void *,
750                 int *numangles,    int **angles,    double **angleforces,
751                 int *numdihedrals, int **dihedrals, double **dihedralforces,
752                 int *numimpropers, int **impropers, double **improperforces,
753                 int *numcterms,    int **cterms, 
754                 int *ctermcols,    int *ctermrows,  double **ctermforces);
755
756   /** 
757    * Write out Angles, Dihedrals, Impropers, and Cross Terms
758    * Forces are in Kcal/mol
759    * (Cross terms pertain to the CHARMM/NAMD CMAP feature, forces are given
760    *  as a 2-D matrix)
761    */
762   int (* write_angles)(void *,
763         int numangles,    const int *angles,    const double *angleforces,
764         int numdihedrals, const int *dihedrals, const double *dihedralforces,
765         int numimpropers, const int *impropers, const double *improperforces,
766         int numcterms,   const int *cterms,    
767         int ctermcols, int ctermrows, const double *ctermforces);
768 #endif
769
770   /**
771    * Retrieve metadata pertaining to QM datasets in this file.
772    */
773   int (* read_qm_metadata)(void *, molfile_qm_metadata_t *metadata);
774
775   /**
776    * Read QM data
777    */
778   int (* read_qm_rundata)(void *, molfile_qm_t *qmdata);
779
780   /**
781    * Read the next timestep from the file.  Return MOLFILE_SUCCESS, or 
782    * MOLFILE_EOF on EOF.  If the molfile_timestep_t or molfile_qm_metadata_t
783    * arguments are NULL, then the coordinate or qm data should be skipped.  
784    * Otherwise, the application must prepare molfile_timestep_t and 
785    * molfile_qm_timestep_t by allocating space for the corresponding 
786    * number of coordinates, orbital wavefunction coefficients, etc.
787    * Since it is common for users to want to load only the final timestep 
788    * data from a QM run, the application may provide any combination of
789    * valid, or NULL pointers for the molfile_timestep_t and 
790    * molfile_qm_timestep_t parameters, depending on what information the
791    * user is interested in.
792    * The natoms and qm metadata parameters exist because some file formats 
793    * cannot determine for themselves how many atoms etc are in a 
794    * timestep; the app must therefore obtain this information elsewhere
795    * and provide it to the plugin.
796    */
797   int (* read_timestep)(void *, int natoms, molfile_timestep_t *,
798                         molfile_qm_metadata_t *, molfile_qm_timestep_t *);
799 #endif
800
801 #if vmdplugin_ABIVERSION > 10
802   int (* read_timestep_metadata)(void *, molfile_timestep_metadata_t *);
803 #endif
804 #if vmdplugin_ABIVERSION > 11
805   int (* read_qm_timestep_metadata)(void *, molfile_qm_timestep_metadata_t *);
806 #endif
807
808 #if vmdplugin_ABIVERSION > 13
809   /**
810    *  Console output, READ-ONLY function pointer.
811    *  Function pointer that plugins can use for printing to the host
812    *  application's text console.  This provides a clean way for plugins
813    *  to send message strings back to the calling application, giving the
814    *  caller the ability to prioritize, buffer, and redirect console messages
815    *  to an appropriate output channel, window, etc.  This enables the use of
816    *  graphical consoles like TkCon without losing console output from plugins.
817    *  If the function pointer is NULL, no console output service is provided
818    *  by the calling application, and the output should default to stdout
819    *  stream.  If the function pointer is non-NULL, all output will be
820    *  subsequently dealt with by the calling application.
821    *
822    *  XXX this should really be put into a separate block of
823    *      application-provided read-only function pointers for any
824    *      application-provided services
825    */
826   int (* cons_fputs)(const int, const char*);
827 #endif
828
829 } molfile_plugin_t;
830
831 #endif
832